Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Channeling.cc
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
27#include "G4Channeling.hh"
28
29#include "Randomize.hh"
30
32#include "G4TouchableHistory.hh"
33#include "G4SystemOfUnits.hh"
34#include "G4LambdacPlus.hh"
36
38G4VDiscreteProcess("channeling"),
39fChannelingID(G4PhysicsModelCatalog::GetModelID("model_channeling")),
40fTimeStepMin(0.),
41fTimeStepMax(0.),
42fTransverseVariationMax(2.E-2 * CLHEP::angstrom),
43k010(G4ThreeVector(0.,1.,0.)),
44fSpin(G4ThreeVector(0.,0.,0.))
45{}
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
53G4ChannelingTrackData* G4Channeling::GetTrackData(const G4Track& aTrack){
54 G4ChannelingTrackData* trackdata =
56 if(trackdata == nullptr){
57 trackdata = new G4ChannelingTrackData();
58 aTrack.SetAuxiliaryTrackInformation(fChannelingID,trackdata);
59 }
60 return trackdata;
61}
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
65void G4Channeling::GetEF(const G4Track& aTrack,
66 G4ThreeVector& pos,
67 G4ThreeVector& out){
68 out = G4ThreeVector((GetMatData(aTrack)->GetEFX()->GetEC(pos)),
69 (GetMatData(aTrack)->GetEFY()->GetEC(pos)),
70 0.);
71}
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
76 G4TouchableHistory* theTouchable = (G4TouchableHistory*)(step->GetTouchable());
77
78 pos -= theTouchable->GetTranslation();
79 pos = ((*theTouchable->GetRotation()).inverse())(pos);
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
84G4bool G4Channeling::UpdateParameters(const G4Track& aTrack){
85
87
88 G4StepPoint* postStepPoint = aTrack.GetStep()->GetPostStepPoint();
89 G4StepPoint* preStepPoint = aTrack.GetStep()->GetPreStepPoint();
90
91 G4ThreeVector posPost = postStepPoint->GetPosition();
92 aLCV->RotateToLattice(posPost);
93 G4ThreeVector posPre = preStepPoint->GetPosition();
94 aLCV->RotateToLattice(posPre);
95
96 G4double integrationLimit = std::fabs(posPost.z() - posPre.z());
97
98 if(integrationLimit > 0.){
99 //----------------------------------------
100 // Check if the crystal is bent
101 //----------------------------------------
102 G4bool isBent = GetMatData(aTrack)->IsBent();
103
104 //----------------------------------------
105 // Get the momentum in the world reference
106 // frame and rotate to the solid reference frame
107 //----------------------------------------
108 G4TouchableHistory* theTouchable = (G4TouchableHistory*)(preStepPoint->GetTouchable());
109 G4ThreeVector momWorld = aTrack.GetStep()->GetPreStepPoint()->GetMomentum();
110 G4ThreeVector mom = (*theTouchable->GetRotation())(momWorld);
111
112 //----------------------------------------
113 // Get the momentum in the solid reference
114 // frame and rotate to the crystal reference frame
115 //----------------------------------------
116 aLCV->RotateToLattice(mom);
117
118 //----------------------------------------
119 // Get the momentum in the crystal reference
120 // frame and rotate to the reference frame
121 // solidal to the bent planes
122 //----------------------------------------
123 if(isBent){
124 PosToLattice(preStepPoint,posPre);
125 G4ThreeVector axis010 = (*theTouchable->GetRotation())(k010);
126 mom.rotate(axis010,-posPre.z()/GetMatData(aTrack)->GetBR(posPre).x());
127 }
128
129 //----------------------------------------
130 // Take the position stored in the track data.
131 // If the particle enters the crystal,
132 // the position in the channel is randomly
133 // generated using a uniform distribution
134 //----------------------------------------
136 if(GetTrackData(aTrack)->GetPosCh().x() == DBL_MAX){
137 G4double posX = G4UniformRand() * GetMatData(aTrack)->GetPot()->GetIntSp(0);
138 G4double posY = G4UniformRand() * GetMatData(aTrack)->GetPot()->GetIntSp(1);
139 pos = G4ThreeVector(posX,posY,0.);
140 }
141 else{
142 pos = GetTrackData(aTrack)->GetPosCh();
143 }
144
145 G4double step=0., stepTot=0.;
146 G4double nud =0., eld =0.;
147 G4double efx =0., efy =0.;
148 G4double nud_temp =0., eld_temp =0.;
149
150 G4double beta = aTrack.GetVelocity()/CLHEP::c_light;
151 G4double Z = GetParticleDefinition(aTrack)->GetPDGCharge();
152
153 const G4double oneSixth = 1./6.;
154 G4ThreeVector posk1,posk2,posk3,posk4,posk5,posk6;
155 G4ThreeVector momk1,momk2,momk3,momk4,momk5,momk6;
156 G4ThreeVector pos_temp, efxy;
157
158 do{
159 //----------------------------------------
160 // Limit the variable step length for the
161 // integration via the selected algorithm
162 // and update variables for the integration
163 //----------------------------------------
164
165 UpdateIntegrationStep(aTrack,mom,step);
166 if(step + stepTot > integrationLimit){
167 step = integrationLimit - stepTot;
168 }
169
170 //----------------------------------------
171 // Function integration algorithm
172 // 4th Order Runge-Kutta
173 //----------------------------------------
174
175 GetEF(aTrack,pos,efxy);
176 posk1 = step / mom.z() * mom;
177 momk1 = step / beta * Z * efxy;
178 if(isBent) momk1.setX(momk1.x() - step * mom.z() * beta / (GetMatData(aTrack)->GetBR(pos)).x());
179
180 GetEF(aTrack,pos_temp = pos + posk1 * 0.5,efxy);
181 posk2 = step / mom.z() * (mom + momk1 * 0.5);
182 momk2 = step / beta * Z * efxy;
183 if(isBent) momk2.setX(momk2.x() - step * mom.z() * beta / (GetMatData(aTrack)->GetBR(pos_temp)).x());
184
185 GetEF(aTrack,pos_temp = pos + posk2 * 0.5,efxy);
186 posk3 = step / mom.z() * (mom + momk2 * 0.5);
187 momk3 = step / beta * Z * efxy;
188 if(isBent) momk3.setX(momk3.x() - step * mom.z() * beta / (GetMatData(aTrack)->GetBR(pos_temp)).x());
189
190 GetEF(aTrack,pos_temp = pos + posk3,efxy);
191 posk4 = step / mom.z() * (mom + momk3);
192 momk4 = step / beta * Z * efxy;
193 if(isBent) momk4.setX(momk4.x() - step * mom.z() * beta / (GetMatData(aTrack)->GetBR(pos_temp)).x());
194
195 pos = pos + oneSixth * (posk1 + 2.*posk2 + 2.*posk3 + posk4);
196 mom = mom + oneSixth * (momk1 + 2.*momk2 + 2.*momk3 + momk4);
197
198 //----------------------------------------
199 // Function integration algorithm
200 // 2th Order Velocity-Verlet
201 //----------------------------------------
202
203 /*
204 GetEF(aTrack,pos,efxy);
205 posk1 = pos + (step * 0.5 / mom.z()) * mom;
206 //momk1 = mom + step * 0.5 / betaZ * efxy;
207 momk1 = mom;
208 if(isBent) momk1.setX(momk1.x() - step * 0.5 * mom.z() * beta / (GetMatData(aTrack)->GetBR(pos)).x());
209
210 GetEF(aTrack,posk1,efxy);
211 pos = pos + (step / momk1.z()) * momk1;
212 //mom = mom + step / betaZ * efxy;
213 mom = mom;
214 if(isBent) mom.setX(mom.x() - step * mom.z() * beta / (GetMatData(aTrack)->GetBR(posk1)).x());
215 */
216
217 //----------------------------------------
218 // Update the total step and the electron
219 // and nuclei density experienced by
220 // the particle during its motion
221 //----------------------------------------
222
223 stepTot += step;
224
225 nud_temp = GetMatData(aTrack)->GetNuD()->GetEC(pos);
226 eld_temp = GetMatData(aTrack)->GetElD()->GetEC(pos);
227
228 if(nud_temp < 0.) {nud_temp = 0.;}
229 if(eld_temp < 0.) {eld_temp = 0.;}
230
231 nud += (step * nud_temp);
232 eld += (step * eld_temp);
233
234 efx += (step * GetMatData(aTrack)->GetEFX()->GetEC(pos));
235 efy += (step * GetMatData(aTrack)->GetEFY()->GetEC(pos));
236
237 } while(stepTot<integrationLimit);
238
239 nud /= stepTot;
240 eld /= stepTot;
241
242 if(nud < 1.E-10) {nud = 1.E-10;}
243 if(eld < 1.E-10) {eld = 1.E-10;}
244
245 GetTrackData(aTrack)->SetNuD(nud);
246 GetTrackData(aTrack)->SetElD(eld);
247
248 GetTrackData(aTrack)->SetEFX(efx);
249 GetTrackData(aTrack)->SetEFY(efy);
250
251 GetTrackData(aTrack)->SetMomCh(mom);
252 GetTrackData(aTrack)->SetPosCh(pos);
253 return true;
254 }
255 else{
256 return false;
257 }
258
259 return false;
260}
261
262//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
263
264G4bool G4Channeling::
265UpdateIntegrationStep(const G4Track& aTrack,
266 G4ThreeVector& mom,
267 G4double& step){
268
269 if(mom.x() != 0.0 || mom.y() != 0.0){
270 double xy2 = mom.x() * mom.x() + mom.y()*mom.y();
271
272 if(xy2!=0.){
273 step = std::fabs(fTransverseVariationMax * GetPre(aTrack)->GetKineticEnergy() / std::pow(xy2,0.5));
274 if(step < fTimeStepMin) step = fTimeStepMin;
275 else{
276 fTimeStepMax = std::sqrt( fTransverseVariationMax * GetPre(aTrack)->GetKineticEnergy()
277 / std::fabs(GetMatData(aTrack)->GetEFX()->GetMax()));
278
279 if(step > fTimeStepMax) step = fTimeStepMax;
280 }
281 }
282 else{
283 step = fTimeStepMin;
284 }
285
286 return true;
287 }
288 else{
289 step = fTimeStepMin;
290 }
291 return false;
292}
293
294//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
295
297GetMeanFreePath(const G4Track& aTrack,
298 G4double, // previousStepSize
300
301 //----------------------------------------
302 // the condition is forced to check if
303 // the volume has a lattice at each step.
304 // if it hasn't, return DBL_MAX
305 //----------------------------------------
306
307 *condition = Forced;
308
309 G4LogicalVolume* aLV = aTrack.GetVolume()->GetLogicalVolume();
311
312 if(G4LogicalCrystalVolume::IsLattice(aLV) == true &&
314 G4double osc_per = GetOscillationPeriod(aTrack);
315 fTimeStepMin = osc_per * 2.E-4;
316 return osc_per * 0.01;
317 }
318 else{
319 GetTrackData(aTrack)->Reset();
320 return DBL_MAX;
321 }
322}
323
324//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
325
327PostStepDoIt(const G4Track& aTrack,
328 const G4Step&){
329
330 //----------------------------------------
331 // check if the volume has a lattice
332 // and if the particle is in channeling.
333 // If it is so, the particle is forced
334 // to follow the channeling plane
335 // direction. If the particle has
336 // dechanneled or exited the crystal,
337 // the outgoing angle is evaluated
338 //----------------------------------------
339
341 G4LogicalVolume* aLV = aTrack.GetVolume()->GetLogicalVolume();
343
344
345 if(G4LogicalCrystalVolume::IsLattice(aLV) == true &&
347
348 G4bool bModifiedTraj = UpdateParameters(aTrack);
349
350 if(bModifiedTraj==true){
351 //----------------------------------------
352 // Get the momentum in the reference frame
353 // solidal to the bent planes and rotate
354 // to the reference frame
355 //----------------------------------------
357 G4ThreeVector momCh = GetTrackData(aTrack)->GetMomCh();
358
359 G4StepPoint* postStepPoint = aTrack.GetStep()->GetPostStepPoint();
360 G4TouchableHistory* theTouchable = (G4TouchableHistory*)(postStepPoint->GetTouchable());
361
362 if(GetMatData(aTrack)->IsBent()){
363 G4ThreeVector posPost = postStepPoint->GetPosition();
364 PosToLattice(postStepPoint,posPost);
365 G4ThreeVector axis010 = (*theTouchable->GetRotation())(k010);
366 momCh.rotate(axis010,posPost.z()/GetMatData(aTrack)->GetBR(posPost).x());
367 }
368
369 //----------------------------------------
370 // Get the momentum in the crystal reference
371 // frame and rotate to the solid reference frame
372 //----------------------------------------
373
374 aLCV->RotateToSolid(momCh);
375
376 //----------------------------------------
377 // Get the momentum in the solid reference
378 // frame and rotate to the world reference frame
379 //----------------------------------------
380 G4ThreeVector mom = ((*theTouchable->GetRotation()).inverse())(momCh);
381
384 }
385 }
386 else{
387 // if the volume has no lattice it resets the density factors
388 GetTrackData(aTrack)->Reset();
389 }
390
391 return &aParticleChange;
392}
393
394//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ Forced
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
const G4int Z[17]
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
Hep3Vector unit() const
double x() const
double y() const
void setX(double)
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:24
G4double GetIntSp(G4int index)
G4double GetEC(G4ThreeVector &)
virtual G4ThreeVector GetBR(G4ThreeVector &v3)
void SetPosCh(G4ThreeVector a3vec)
void SetNuD(G4double aDouble)
void SetElD(G4double aDouble)
void SetEFY(G4double aDouble)
void SetEFX(G4double aDouble)
void SetMomCh(G4ThreeVector a3vec)
void PosToLattice(G4StepPoint *step, G4ThreeVector &)
Definition: G4Channeling.cc:75
virtual ~G4Channeling()
Definition: G4Channeling.cc:49
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double GetOscillationPeriod(const G4Track &aTrack)
Definition: G4Channeling.hh:87
const G4ThreeVector & RotateToSolid(G4ThreeVector &dir) const
static G4bool IsLattice(G4LogicalVolume *aLV)
const G4ThreeVector & RotateToLattice(G4ThreeVector &dir) const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
void Initialize(const G4Track &) override
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetPDGCharge() const
const G4VTouchable * GetTouchable() const
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetPosition() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
const G4RotationMatrix * GetRotation(G4int depth=0) const
const G4ThreeVector & GetTranslation(G4int depth=0) const
G4double GetVelocity() const
G4VPhysicalVolume * GetVolume() const
void SetAuxiliaryTrackInformation(G4int id, G4VAuxiliaryTrackInformation *info) const
Definition: G4Track.cc:209
G4VPhysicalVolume * GetNextVolume() const
G4VAuxiliaryTrackInformation * GetAuxiliaryTrackInformation(G4int id) const
Definition: G4Track.cc:229
const G4Step * GetStep() const
G4LogicalVolume * GetLogicalVolume() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
Definition: DoubConv.h:17
#define DBL_MAX
Definition: templates.hh:62