Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VChannelingFastSimCrystalData.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Author: Alexei Sytov
27// Co-author: Gianfranco PaternĂ² (modifications & testing)
28// On the base of the CRYSTALRAD realization of scattering model:
29// A. I. Sytov, V. V. Tikhomirov, and L. Bandiera PRAB 22, 064601 (2019)
30
31#ifndef G4VChannelingFastSimCrystalData_h
32#define G4VChannelingFastSimCrystalData_h 1
33
34#include "globals.hh"
35#include "G4ios.hh"
36#include "G4ThreeVector.hh"
37#include "Randomize.hh"
38#include "G4LogicalVolume.hh"
39#include "G4Material.hh"
40#include "G4VSolid.hh"
41#include <unordered_map>
43
45
46/** \file G4VChannelingFastSimCrystalData.hh
47* \brief Definition of the G4VChannelingFastSimCrystalData class
48* The class contains the data and properties related to the crystal lattice as well as
49* functions to simulate of important physical processes, i.e. coulomb scattering on
50* screened atomic potential, on single electrons and ionization energy losses;
51* functions of electric fields, nuclear and electron densities and minimum energy
52* of ionization (the corresponding interpolation coefficients are in
53* G4ChannelingFastSimInterpolation).
54* The functions related to the crystal geometry (transformation of coordinates and angles
55* from the reference system of the bounding box of the local volume to
56* the crystal lattice co-rotating reference system and vice versa) and
57* initialization function SetMaterialProperties are created as virtual to make
58* material data input and geometry functions flexible for modification.
59*/
60
62public:
63
66
67 ///electric fields produced by crystal lattice
68 G4double Ex(G4double x,G4double y) {return (fElectricFieldX->GetIF(x,y))*(-fZ2/fPV);}
69 G4double Ey(G4double x,G4double y) {return (fElectricFieldY->GetIF(x,y))*(-fZ2/fPV);}
70
71 ///electron density function
73 {
74 G4double nel0=fElectronDensity->GetIF(x,y);
75 if(nel0<0.) {nel0=0.;}//exception, errors of interpolation functions
76 return nel0;
77 }
78 ///minimum energy of ionization function
81 ///nuclear density function (normalized to average nuclear density)
83 {return std::abs(fNucleiDensity[ielement]->GetIF(x,y));}
84 //abs to describe exception, errors of interpolation functions,
85 //don't put it =0, otherwise division on 0 in CoulombAtomicScattering
86
87 ///Calculate the value of the Lindhard angle (!!! the value for a straight crystal)
89 ///Calculate the value of the Lindhard angle (!!! the value for a straight crystal)
90 G4double GetLindhardAngle();//return the Lindhard angle value calculated in
91 //SetParticleProperties
92
93 ///Calculate simulation step (standard value for channeling particles and
94 ///reduced value for overbarrier particles)
96 ///Calculate maximal simulation step (standard value for channeling particles)
98
99 ///get particle velocity/c
100 G4double GetBeta(){return fBeta;}
101
103 G4int GetModel() {return iModel;}//=1 for planes, =2 for axes
104
105 ///get bending angle of the crystal planes/axes
106 ///(default BendingAngle=0 => straight crystal);
108 ///fBendingAngle MAY BE NOT THE SAME AS THE BENDING ANGLE OF THE CRYSTAL VOLUME:
109 ///THE VOLUME OF A BENT CRYSTAL MAY BE G4Box, while the planes/axes inside may be bent
110
112
113 ///get crystal curvature
114 ///for crystalline undulator the curvature is a function, otherwise it's a constant
116 {return fCU ? //select between a crystalline undulator (CU) and a bent crystal
117 ( fImportCrystalGeometry ? //select between a realistic and an ideal CU
118 fVecCUCurv[fCUID].Value(z) :
119 -fCUK2*GetCUx(z)) :
120 fCurv;}
121
122 ///get crystalline undulator wave function
124 {return fImportCrystalGeometry ? //select between a realistic and an ideal CU
125 fVecCUx[fCUID].Value(z) :
126 fCUAmplitude*std::cos(fCUK*z+fCUPhase);}
127
128 ///get crystalline undulator wave 1st derivative function
130 {return fCU ? //select between a crystalline undulator (CU) and a bent crystal
131 ( fImportCrystalGeometry ? //select between a realistic and an ideal CU
132 fVecCUtetax[fCUID].Value(z) :
133 -fCUAmplitudeK*std::sin(fCUK*z+fCUPhase)) :
134 0.;}
135
136 ///find and upload crystal lattice input files, calculate all the basic values
137 ///(to do only once)
138 virtual void SetMaterialProperties(const G4Material* crystal,
139 const G4String &lattice,
140 const G4String &filePath) = 0;
141
142 ///set geometry parameters from current logical volume
143 void SetGeometryParameters(const G4LogicalVolume *crystallogic);
144
145 ///set bending angle of the crystal planes/axes
146 ///(default fBendingAngle=0 => straight crystal);
147 ///only non-negative values! crystal is bent in the positive direction of x
148 void SetBendingAngle(G4double tetab, const G4LogicalVolume *crystallogic);
149 ///fBendingAngle MAY BE NOT THE SAME AS THE BENDING ANGLE OF THE CRYSTAL VOLUME
150 ///THE VOLUME OF A BENT CRYSTAL MAY BE G4Box, while the planes/axes inside may be bent
151
152 ///set miscut angle (default fMiscutAngle=0), acceptable range +-1 mrad,
153 ///otherwise geometry routines may be unstable
154 void SetMiscutAngle(G4double tetam, const G4LogicalVolume *crystallogic);
155
156 ///set crystalline undulator parameters: amplitude, period and phase
157 /// (default: all 3 value = 0)
158 /// function to use in Detector Construction
160 G4double period,
161 G4double phase,
162 const G4LogicalVolume *crystallogic);
163
164 ///set importing geometry of a crystalline undulator from file for specific volume
165 void SetCrystallineUndulatorParameters(const G4LogicalVolume *crystallogic,
166 const G4String &filename = "CUgeometry.dat");
167
168 ///set crystalline undulator parameters (internal function of the model)
169 ///for convenience we put amplitude, period and phase in a G4ThreeVector
170 void SetCUParameters(const G4ThreeVector &amplitudePeriodPhase,
171 const G4LogicalVolume *crystallogic);
172
173 ///recalculate all the important values
174 ///(to do both at the trajectory start and after energy loss)
176 G4double mp,
177 G4double charge,
178 const G4String& particleName);
179
180 ///calculate the coordinates in the co-rotating reference system
181 ///within a channel (periodic cell)
182 ///(connected with crystal planes/axes either bent or straight)
184
185 ///calculate the coordinates in the Box reference system
186 ///(connected with the bounding box of the volume)
188
189 ///change the channel if necessary, recalculate x o y
191
192 ///return correction of the longitudinal coordinate
193 /// (along current plane/axis vs "central plane/axis")
195
196 ///calculate the horizontal angle in the co-rotating reference system
197 ///within a channel (periodic cell)
198 ///(connected with crystal planes/axes either bent or straight)
200
201 ///calculate the horizontal angle in the Box reference system
202 ///(connected with the bounding box of the volume)
204
205 ///auxialiary function to transform the horizontal angle
207
208 ///multiple and single scattering on screened potential
210 G4double effectiveStep,
211 G4double step,
212 G4int ielement);
213 ///multiple and single scattering on electrons
215 G4double electronDensity,
216 G4double step);
217 ///ionization losses
219
220 void SetVerbosity(G4int ver){fVerbosity = ver;}
221
222protected:
223 ///classes containing interpolation coefficients
224 //horizontal electric field data
226 //vertical electric field data
228 //electron density data
230 //minimal energy of ionization data
232 //nuclear density distributions data
233 std::vector <G4ChannelingFastSimInterpolation*> fNucleiDensity;
234
235 ///values related to the crystal geometry
236 G4ThreeVector fHalfDimBoundingBox;//bounding box half dimensions
237 G4int fBent=0;//flag of bent crystal,
238 //=0 for straight and =1 for bent, by default straight crystal
239
240 G4double fBendingAngle=0.;// angle of bending of the crystal planes/axes
241 //inside the crystal volume
242 //MAY BE NOT THE SAME AS THE BENDING ANGLE OF THE CRYSTAL VOLUME
243 //THE VOLUME OF A BENT CRYSTAL MAY BE G4Box,
244 //while the planes/axes inside may be bent
245 G4double fBendingR = 0.; // bending radius of the crystal planes/axes
246 G4double fBending2R=0.; // =2*fBendingR
247 G4double fBendingRsquare=0.; // =fBendingR**2
248 G4double fCurv=0.; //=1/fBendingR bending curvature of the crystal planes/axes
249
250 G4double fMiscutAngle = 0.;// miscut angle, can be of either sign or 0;
251 //safe values |ThetaMiscut|<0.001
252 G4double fCosMiscutAngle=1.;// = std::cos(fMiscutAngle), to economy operations
253 G4double fSinMiscutAngle=0.;// = std::sin(fMiscutAngle), to economy operations
254
255 G4double fCorrectionZ = 1.;//correction of the longitudinal coordinate
256 //(along current plane/axis vs "central plane/axis"), 1 is default value
257 //(for "central plane/axis" or a straight crystal)
258
259 G4bool fCU = false;//flag of crystalline undulator geometry
260 //(periodically bent crystal)
261 G4double fCUAmplitude=0.; //Amplitude of a crystalline undulator
262 G4double fCUK=0.; //2*pi/period of a crystalline undulator
263 G4double fCUPhase=0.;//Phase of a crystalline undulator
264 G4double fCUAmplitudeK=0.;//fCUAmplitude*fCUK
265 G4double fCUK2=0.; //fCUK^2
266
267 ///values related to the crystal lattice
268 G4int fNelements=1;//number of nuclear elements in a crystal
269 G4int iModel=1;// model type (iModel=1 for interplanar potential,
270 //iModel=2 for the interaxial one)
271
272 G4double fVmax=0; // the height of the potential well
273 G4double fVmax2=0; // =2*fVmax
274 G4double fVMinCrystal=0;// non-zero minimal potential inside the crystal,
275 // necessary for angle recalculation for entrance/exit
276 //through the crystal lateral surface
277
278 G4double fChangeStep=0;// fChannelingStep = fChangeStep/fTetaL
279
280 std::vector <G4double> fI0; //Mean excitation energy
281
282 std::vector <G4double> fRF;//Thomas-Fermi screening radius
283
284 ///angles necessary for multiple and single coulomb scattering
285
286 //minimal scattering angle by coulomb scattering on nuclei
287 //defined by shielding by electrons
288 std::vector <G4double> fTeta10;//(in the Channeling model
289 //teta1=fTeta10/fPz*(1.13+fK40/vz**2)
290 //maximal scattering angle by coulomb scattering on nuclei defined by nucleus radius
291 std::vector <G4double> fTetamax0;//(in the Channeling model tetamax=fTetamax0/fPz)
292 std::vector <G4double> fTetamax2;//=tetamax*tetamax
293 std::vector <G4double> fTetamax12;//=teta1*teta1+tetamax*tetamax
294 std::vector <G4double> fTeta12; //= teta1*teta1
295
296 ///coefficients necessary for multiple and single coulomb scattering
297 std::vector <G4double> fK20; //a useful coefficient, fK2=fK20/fPV/fPV
298 std::vector <G4double> fK2; //a useful coefficient,
299 //fK2=(fZ2*alpha*hdc)**2*4.*pi*fN0*(fZ1/fPV)**2
300 std::vector <G4double> fK40; //a useful coefficient, fK40=3.76D0*(alpha*fZ1)**2
301 G4double fK30=0;//a useful coefficient, fK3=fK30/fPV/fPV
302 G4double fK3=0;//a useful coefficient, fK3=2.*pi*alpha*hdc/electron_mass_c2/(fPV)**2
303
304 std::vector <G4double> fKD; //a useful coefficient for dE/dx
305 std::vector <G4double> fLogPlasmaEdI0; //item of delta-correction of ionization loss
306
307 ///coefficients for multiple scattering suppression
308 std::vector <G4double> fPu11;//a useful coefficient for exponent containing u1
309 std::vector <G4double> fPzu11;//a useful coefficient for exponent containing u1
310 std::vector <G4double> fBB;//a useful coefficient
311 std::vector <G4double> fE1XBbb;//a useful coefficient
312 std::vector <G4double> fBBDEXP;//a useful coefficient
313
314 //Variable to control printout
316
317private:
318
319 ///variables for realistic crystalline undulator
320 // flag of custom crystal geometry uploaded from a file
321 G4bool fImportCrystalGeometry = false;
322
323 //crystalline undulator wave function data (realistic crystalline undulator)
324 std::vector <G4PhysicsLinearVector> fVecCUx;
325 //crystalline undulator 1st derivative wave function data (realistic crystalline undulator)
326 std::vector <G4PhysicsLinearVector> fVecCUtetax;
327 //crystalline undulator curvature (realistic crystalline undulator)
328 std::vector <G4PhysicsLinearVector> fVecCUCurv;
329
330 std::unordered_map<G4int, G4ThreeVector> fMapCUAmplitudePeriodPhase;//the map of
331 //AmplitudePeriodPhase
332 //for different logical volumes
333
334 std::unordered_map<G4int, G4bool> fMapImportCrystalGeometry;//the map of
335 //fImportCrystalGeometry flag
336 //for different logical volumes
337
338 std::unordered_map<G4int, G4int> fMapCUID;//the map of the crystalline undulator ID
339
340 G4int fCUID = 0; //current logical volume ID for the crystalline undulator
341
342 //exponential integral
343 G4double expint(G4double x);
344
345 ///private variables
346
347 std::unordered_map<G4int, G4double> fMapBendingAngle;//the map fBendingAngle
348 //for different logical volumes
349
350 std::unordered_map<G4int, G4double> fMapMiscutAngle;//the map fMiscutAngle
351 //for different logical volumes
352
353
354 G4double fChannelingStep=0;// simulation step under the channeling conditions =
355 //channeling oscillation length/fNsteps
356 // channeling oscillation length: Biryukov book Eq. (1.24)
357
358 ///energy depended values
359 G4double fPz=0; // particle momentum absolute value
360 G4double fPV=0; // pv
361 G4double fTetaL=0; //Lindhard angle
362 G4double fBeta=0; //particle (velocity/c)
363 G4double fV2=0; // particle (velocity/c)^2
364 G4double fGamma=0; //Lorentz factor
365 G4double fMe2Gamma=0; // me^2*fGamma
366 G4double fTmax=0; // max ionization losses
367
368 ///particle properties flags
369 G4String fParticleName = "";
370 G4double fZ2=0; //particle charge
371
372};
373
374#endif
375
Definition of the G4ChannelingFastSimInterpolation class The class includes spline interpolation coef...
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
virtual G4double AngleXShift(G4double z)=0
auxialiary function to transform the horizontal angle
G4double IonizationLosses(G4double dz, G4int ielement)
ionization losses
std::vector< G4double > fK20
coefficients necessary for multiple and single coulomb scattering
G4double GetLindhardAngle()
Calculate the value of the Lindhard angle (!!! the value for a straight crystal).
G4double GetBeta()
get particle velocity/c
std::vector< G4ChannelingFastSimInterpolation * > fNucleiDensity
void SetGeometryParameters(const G4LogicalVolume *crystallogic)
set geometry parameters from current logical volume
void SetCUParameters(const G4ThreeVector &amplitudePeriodPhase, const G4LogicalVolume *crystallogic)
G4double GetMaxSimulationStep(G4double etotal, G4double mass, G4double charge)
Calculate maximal simulation step (standard value for channeling particles).
G4ChannelingFastSimInterpolation * fMinIonizationEnergy
G4ChannelingFastSimInterpolation * fElectricFieldY
virtual G4ThreeVector ChannelChange(G4double &x, G4double &y, G4double &z)=0
change the channel if necessary, recalculate x o y
virtual G4ThreeVector CoordinatesFromBoxToLattice(const G4ThreeVector &pos0)=0
virtual G4double AngleXFromBoxToLattice(G4double tx, G4double z)=0
G4double ElectronDensity(G4double x, G4double y)
electron density function
void SetCrystallineUndulatorParameters(G4double amplitude, G4double period, G4double phase, const G4LogicalVolume *crystallogic)
G4ThreeVector fHalfDimBoundingBox
values related to the crystal geometry
G4double GetSimulationStep(G4double tx, G4double ty)
void SetBendingAngle(G4double tetab, const G4LogicalVolume *crystallogic)
virtual void SetMaterialProperties(const G4Material *crystal, const G4String &lattice, const G4String &filePath)=0
virtual G4double AngleXFromLatticeToBox(G4double tx, G4double z)=0
G4double MinIonizationEnergy(G4double x, G4double y)
minimum energy of ionization function
G4ChannelingFastSimInterpolation * fElectricFieldX
classes containing interpolation coefficients
G4double Ex(G4double x, G4double y)
electric fields produced by crystal lattice
G4ChannelingFastSimInterpolation * fElectronDensity
G4int fNelements
values related to the crystal lattice
G4ThreeVector CoulombElectronScattering(G4double eMinIonization, G4double electronDensity, G4double step)
multiple and single scattering on electrons
std::vector< G4double > fTeta10
angles necessary for multiple and single coulomb scattering
G4double NuclearDensity(G4double x, G4double y, G4int ielement)
nuclear density function (normalized to average nuclear density)
G4double GetCUx(G4double z)
get crystalline undulator wave function
virtual G4ThreeVector CoordinatesFromLatticeToBox(const G4ThreeVector &pos)=0
std::vector< G4double > fPu11
coefficients for multiple scattering suppression
G4ThreeVector CoulombAtomicScattering(G4double effectiveStep, G4double step, G4int ielement)
multiple and single scattering on screened potential
G4double GetCUtetax(G4double z)
get crystalline undulator wave 1st derivative function
void SetMiscutAngle(G4double tetam, const G4LogicalVolume *crystallogic)
void SetParticleProperties(G4double etotal, G4double mp, G4double charge, const G4String &particleName)