Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChannelingFastSimModel Class Reference

#include <G4ChannelingFastSimModel.hh>

+ Inheritance diagram for G4ChannelingFastSimModel:

Public Member Functions

 G4ChannelingFastSimModel (const G4String &, G4Region *)
 
 G4ChannelingFastSimModel (const G4String &)
 
 ~G4ChannelingFastSimModel ()
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 – IsApplicable
 
G4bool ModelTrigger (const G4FastTrack &) override
 – ModelTrigger
 
void DoIt (const G4FastTrack &, G4FastStep &) override
 – User method DoIt
 
void Input (const G4Material *crystal, const G4String &lattice)
 special functions
 
void RadiationModelActivate ()
 
G4ChannelingFastSimCrystalDataGetCrystalData ()
 
G4BaierKatkovGetRadiationModel ()
 
G4bool GetIfRadiationModelActive ()
 
void SetLowKineticEnergyLimit (G4double ekinetic, const G4String &particleName)
 set cuts
 
void SetLindhardAngleNumberHighLimit (G4double angleNumber, const G4String &particleName)
 
void SetDefaultLowKineticEnergyLimit (G4double ekinetic)
 
void SetDefaultLindhardAngleNumberHighLimit (G4double angleNumber)
 
void SetMaxPhotonsProducedPerStep (G4double nPhotons)
 
G4double GetLowKineticEnergyLimit (const G4String &particleName)
 get cuts
 
G4double GetLindhardAngleNumberHighLimit (const G4String &particleName)
 
G4double GetLowKineticEnergyLimit (G4int particleDefinitionID)
 
G4double GetLindhardAngleNumberHighLimit (G4int particleDefinitionID)
 
G4int GetMaxPhotonsProducedPerStep ()
 get the maximal number of photons that can be produced per fastStep
 
- Public Member Functions inherited from G4VFastSimulationModel
 G4VFastSimulationModel (const G4String &aName)
 
 G4VFastSimulationModel (const G4String &aName, G4Envelope *, G4bool IsUnique=FALSE)
 
virtual ~G4VFastSimulationModel ()=default
 
virtual G4bool AtRestModelTrigger (const G4FastTrack &)
 
virtual void AtRestDoIt (const G4FastTrack &, G4FastStep &)
 
virtual void Flush ()
 
const G4String GetName () const
 
G4bool operator== (const G4VFastSimulationModel &) const
 

Detailed Description

Definition at line 50 of file G4ChannelingFastSimModel.hh.

Constructor & Destructor Documentation

◆ G4ChannelingFastSimModel() [1/2]

G4ChannelingFastSimModel::G4ChannelingFastSimModel ( const G4String & modelName,
G4Region * envelope )

Definition at line 40 of file G4ChannelingFastSimModel.cc.

42: G4VFastSimulationModel(modelName, envelope)
43{
44}
G4VFastSimulationModel(const G4String &aName)

◆ G4ChannelingFastSimModel() [2/2]

G4ChannelingFastSimModel::G4ChannelingFastSimModel ( const G4String & modelName)

Definition at line 48 of file G4ChannelingFastSimModel.cc.

49: G4VFastSimulationModel(modelName)
50{
51
52}

◆ ~G4ChannelingFastSimModel()

G4ChannelingFastSimModel::~G4ChannelingFastSimModel ( )

Definition at line 56 of file G4ChannelingFastSimModel.cc.

57{
58}

Member Function Documentation

◆ DoIt()

void G4ChannelingFastSimModel::DoIt ( const G4FastTrack & fastTrack,
G4FastStep & fastStep )
overridevirtual

– User method DoIt

Implements G4VFastSimulationModel.

Definition at line 129 of file G4ChannelingFastSimModel.cc.

131{
132 G4double etotal;//particle total energy
133 G4double etotalPreStep;//etotal at the previous step
134 G4double etotalToSetParticleProperties;//etotal value at which
135 //SetParticleProperties is called
136 G4double mass; //particle mass
137 G4double charge;//particle charge
138 G4double tGlobal; //global time
139 G4double tGlobalPreStep; //global time at the previous step
140 G4ThreeVector xyz0;// the coordinates in the local reference system of the volume
141 G4ThreeVector xyz0PreStep;// xyz at the previous step
142 G4ThreeVector xyz;// the coordinates in the co-rotating reference system within
143 //a channel (elementary periodic cell)
144 G4double x,y,z; // the coordinates in the co-rotating reference system within
145 //a channel (elementary periodic cell)
146 G4double tx0,ty0; // the angles in the local reference system of the volume
147 G4double tx,ty; // the angles in the co-rotating reference system within
148 //a channel (elementary periodic cell)
149 G4double txPreStep,tyPreStep;// tx,ty at the previous step
150 G4ThreeVector momentumDirection;
151 G4ThreeVector scatteringAnglesAndEnergyLoss;//output of scattering functions
152 G4double lindhardAngleNumberHighLimit0; //current high limit of the angle expressed in
153 //[Lindhard angle] units
154
155 //coordinates in Runge-Kutta calculations
156 G4double x1=0.,x2=0.,x3=0.,x4=0.,y1=0.,y2=0.,y3=0.,y4=0.;
157 //angles in Runge-Kutta calculations
158 G4double tx1=0.,tx2=0.,tx3=0.,tx4=0.,ty1=0.,ty2=0.,ty3=0.,ty4=0.;
159 //variables in Runge-Kutta calculations
160 G4double kvx1=0.,kvx2=0.,kvx3=0.,kvx4=0.,kvy1=0.,kvy2=0.,kvy3=0.,kvy4=0.;
161 //simulation step along z (internal step of the model) and its parts
162 G4double dz,dzd3,dzd8;//dzd3 = dz/3; dzd8 = dz/8;
163 //simulation step along the momentum direction
164 G4double momentumDirectionStep;
165 //effective simulation step (taking into account nuclear density along the trajectory)
166 G4double effectiveStep;
167
168 // flag, if Inside(xyz0) switches to kInside
169 G4bool inside = false;
170
171 G4LogicalVolume* crystallogic = fastTrack.GetEnvelopeLogicalVolume();
172 fCrystalData->SetGeometryParameters(crystallogic);
173
174 //set the max number of secondaries (photons) that can be added at this fastStep
175 if (fRad)
176 {
177 fastStep.SetNumberOfSecondaryTracks(fMaxPhotonsProducedPerStep);
178 //reseting the BaierKatkov integral to start it with the new trajectory
179 fBaierKatkov->ResetRadIntegral();//to avoid any memory from the previous trajectory
180 }
181
182 etotal = fastTrack.GetPrimaryTrack()->GetTotalEnergy();
183 mass = fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetPDGMass();
184 charge = fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetPDGCharge();
185
186 // we need to distunguish only charge particles, either leptons or hadrons
187 G4bool hadron =
189
190 lindhardAngleNumberHighLimit0 =
192 GetParticleDefinition()->GetParticleDefinitionID());
193
194 //set fCrystalData parameters depending on the particle parameters
195 fCrystalData->SetParticleProperties(etotal, mass, charge, hadron);
196
197 //global time
198 tGlobal = fastTrack.GetPrimaryTrack()->GetGlobalTime();
199
200 //coordinates in the co-rotating reference system within a channel
201 xyz0= fastTrack.GetPrimaryTrackLocalPosition();
202 xyz = fCrystalData->CoordinatesFromBoxToLattice(xyz0);
203 x=xyz.x();
204 y=xyz.y();
205 z=xyz.z();
206
207 momentumDirection=fastTrack.GetPrimaryTrackLocalDirection();
208 //angle in the co-rotating reference system within a channel
209 //(!!! ONLY FORWARD DIRECTION, momentumDirection.getZ()>0,
210 //valid for high energies defined by the standard energy cuts)
211 tx0 = std::atan(momentumDirection.x()/momentumDirection.z());
212 ty0 = std::atan(momentumDirection.y()/momentumDirection.z());
213
214 //angles in the co-rotating reference system within a channel
215 tx = fCrystalData->AngleXFromBoxToLattice(tx0,z);
216 ty = ty0;
217
218 etotalToSetParticleProperties = etotal*0.999;
219 G4bool inCrystal=true;//flag necessary to escape the cycle (at inCrystal=0;)
220 //do calculations until the particle is inside the volume
221 do
222 {
223 //remember the global time before the next step dz
224 tGlobalPreStep=tGlobal;
225 //remember the coordinates before the next step dz
226 xyz0PreStep = xyz0;
227 //remember the angles and the total energy before the step dz
228 txPreStep = tx;
229 tyPreStep = ty;
230 etotalPreStep = etotal;
231
232 dz = fCrystalData->GetSimulationStep(tx,ty);
233 dzd3=dz/3;
234 dzd8=dz/8;
235
236 //trajectory calculation:
237 //Runge-Cutt "3/8"
238 //fCrystalData->GetCurv(z)*fCrystalData->GetCorrectionZ() is due to dependence of
239 //the radius on x; GetCurv gets 1/R for the central ("central plane/axis")
240
241 //first step
242 kvx1=fCrystalData->Ex(x,y);
243 x1=x+tx*dzd3;
244 tx1=tx+(kvx1-fCrystalData->GetCurv(z)*fCrystalData->GetCorrectionZ())*dzd3;
245 if (fCrystalData->GetModel()==2)
246 {
247 kvy1=fCrystalData->Ey(x,y);
248 y1=y+ty*dzd3;
249 ty1=ty+kvy1*dzd3;
250 }
251
252 //second step
253 kvx2=fCrystalData->Ex(x1,y1);
254 x2=x-tx*dzd3+tx1*dz;
255 tx2=tx-(kvx1-fCrystalData->GetCurv(z)*fCrystalData->GetCorrectionZ())*dzd3+
256 (kvx2-fCrystalData->GetCurv(z)*fCrystalData->GetCorrectionZ())*dz;
257 if (fCrystalData->GetModel()==2)
258 {
259 kvy2=fCrystalData->Ey(x1,y1);
260 y2=y-ty*dzd3+ty1*dz;
261 ty2=ty-kvy1*dzd3+kvy2*dz;
262 }
263
264 //third step
265 kvx3=fCrystalData->Ex(x2,y2);
266 x3=x+(tx-tx1+tx2)*dz;
267 tx3=tx+(kvx1-kvx2+kvx3-fCrystalData->GetCurv(z)*fCrystalData->GetCorrectionZ())*dz;
268 if (fCrystalData->GetModel()==2)
269 {
270 kvy3=fCrystalData->Ey(x2,y2);
271 y3=y+(ty-ty1+ty2)*dz;
272 ty3=ty+(kvy1-kvy2+kvy3)*dz;
273 }
274
275 //fourth step
276 kvx4=fCrystalData->Ex(x3,y3);
277 x4=x+(tx+3.*tx1+3.*tx2+tx3)*dzd8;
278 tx4=tx+(kvx1+3.*kvx2+3.*kvx3+kvx4)*dzd8-
279 fCrystalData->GetCurv(z)*fCrystalData->GetCorrectionZ()*dz;
280 if (fCrystalData->GetModel()==2)
281 {
282 kvy4=fCrystalData->Ey(x3,y3);
283 y4=y+(ty+3.*ty1+3.*ty2+ty3)*dzd8;
284 ty4=ty+(kvy1+3.*kvy2+3.*kvy3+kvy4)*dzd8;
285 }
286 else
287 {
288 y4 =y+ty*dz;
289 ty4=ty;
290 }
291
292 x=x4;
293 tx=tx4;
294 y=y4;
295 ty=ty4;
296
297 z+=dz*fCrystalData->GetCorrectionZ();//motion along the z coordinate
298 //("central plane/axis", no current plane/axis)
299
300 xyz = fCrystalData->ChannelChange(x,y,z);
301 x=xyz.x();
302 y=xyz.y();
303 z=xyz.z();
304
305 //the coordinates in the local reference system of the volume
306 //this vector will be used in the cycle escape condition and
307 //in the radiation model (if activated)
308 xyz0=fCrystalData->CoordinatesFromLatticeToBox(xyz);
309
310 momentumDirectionStep=
311 dz*std::sqrt(1+std::pow(std::tan(tx),2)+std::pow(std::tan(ty),2));
312 tGlobal+=momentumDirectionStep/(fCrystalData->GetBeta())/CLHEP::c_light;
313
314 //default scattering and energy loss 0
315 scatteringAnglesAndEnergyLoss = G4ThreeVector(0.,0.,0.);
316
317 //calculate separately for each element of the crystal
318 for (G4int i = 0; i < fCrystalData->GetNelements(); i++)
319 {
320 //effective step taking into account nuclear density along the trajectory
321 effectiveStep = momentumDirectionStep*fCrystalData->NuclearDensity(x,y,i);
322 //Coulomb scattering on screened atomic potential (both multiple and single)
323 scatteringAnglesAndEnergyLoss += fCrystalData->
324 CoulombAtomicScattering(effectiveStep,momentumDirectionStep,i);
325
326 //Amorphous part of ionization energy losses
327 etotal-=fCrystalData->IonizationLosses(momentumDirectionStep, i);
328 }
329 //electron scattering and coherent part of ionization energy losses
330 scatteringAnglesAndEnergyLoss += fCrystalData->CoulombElectronScattering(
331 fCrystalData->MinIonizationEnergy(x,y),
332 fCrystalData->ElectronDensity(x,y),
333 momentumDirectionStep);
334 tx += scatteringAnglesAndEnergyLoss.x();
335 ty += scatteringAnglesAndEnergyLoss.y();
336 etotal -= scatteringAnglesAndEnergyLoss.z();
337
338 // recalculate the energy depended parameters
339 //(only if the energy decreased enough, not at each step)
340 if (etotalToSetParticleProperties>etotal)
341 {
342 fCrystalData->SetParticleProperties(etotal, mass, charge, hadron);
343 etotalToSetParticleProperties = etotal*0.999;
344 }
345
346 //chain of conditions to escape the cycle
347 // if Inside(xyz0)==kInside has been already true
348 //(a particle has been inside the crystal)
349 if (inside)
350 {
351 // if low energy
352 if (etotal-mass<=GetLowKineticEnergyLimit(fastTrack.GetPrimaryTrack()->
353 GetParticleDefinition()->
354 GetParticleDefinitionID()))
355 {inCrystal = false;}//escape the cycle
356 //check if the angle w.r.t. the axes or planes is too high =>
357 //return to standard Geant4:
358 else if (fCrystalData->GetModel()==1) //1D model, field of planes
359 {
360 //if the angle w.r.t. the planes is too high
361 if (std::abs(tx) >=
362 lindhardAngleNumberHighLimit0*fCrystalData->GetLindhardAngle())
363 {inCrystal = false;}//escape the cycle
364 }
365 else if (fCrystalData->GetModel()==2) //2D model, field of axes
366 {
367 //if the angle w.r.t. the axes is too high
368 if (std::sqrt(tx*tx+ty*ty) >= lindhardAngleNumberHighLimit0*
369 fCrystalData->GetLindhardAngle())
370 {inCrystal = false;}//escape the cycle
371 }
372
373 //radiation production & radiation energy losses
374 //works only if the radiation model is activated
375 if (fRad)
376 {
377 //back to the local reference system of the volume
378 tx0 = fCrystalData->AngleXFromLatticeToBox(tx,z);
379 ty0 = ty;
380 //xyz0 was calculated above
381
382 //running the radiation model and checking if a photon has been emitted
383 if(fBaierKatkov->DoRadiation(etotal,mass,
384 tx0,ty0,
385 scatteringAnglesAndEnergyLoss.x(),
386 scatteringAnglesAndEnergyLoss.y(),
387 momentumDirectionStep,tGlobal,xyz0,
388 crystallogic->
389 GetSolid()->
390 Inside(xyz0)!=kInside&&inCrystal))
391 // also it was checked if the particle is escaping the volume
392 // calculate the radiation integral immidiately in this case
393 {
394 //a photon has been emitted!
395 //shift the particle back into the radiation point
396 etotal = fBaierKatkov->GetParticleNewTotalEnergy();
397 tx0 = fBaierKatkov->GetParticleNewAngleX();
398 ty0 = fBaierKatkov->GetParticleNewAngleY();
399 tGlobal = fBaierKatkov->GetNewGlobalTime();
400 xyz0 = fBaierKatkov->GetParticleNewCoordinateXYZ();
401
402 //add secondary photon
403 fBaierKatkov->GeneratePhoton(fastStep);
404
405 //particle energy was changed
406 fCrystalData->SetParticleProperties(etotal, mass, charge, hadron);
407
408 //coordinates in the co-rotating reference system within a channel
409 xyz = fCrystalData->CoordinatesFromBoxToLattice(xyz0);
410 x=xyz.x();
411 y=xyz.y();
412 z=xyz.z();
413
414 //angles in the co-rotating reference system within a channel
415 tx = fCrystalData->AngleXFromBoxToLattice(tx0,z);
416 ty = ty0;
417 }
418 }
419
420 //precise check if the particle is escaping the volume
421 if (crystallogic->GetSolid()->
422 Inside(xyz0)!=kInside)
423 {
424 //one step back to remain inside the volume
425 //after the escape of the volume
426 tGlobal = tGlobalPreStep;
427 xyz0 = xyz0PreStep;
428 tx = txPreStep;
429 ty = tyPreStep;
430 etotal = etotalPreStep;
431 z-=dz*fCrystalData->GetCorrectionZ();
432 // change the flag => this particle will not enter
433 // the model before escape this volume
434
435 inCrystal = false; //escape the cycle
436 }
437 }
438 else
439 {
440 // if Inside(xyz0)==kInside we can enable checking of particle escape
441 if (crystallogic->GetSolid()->
442 Inside(xyz0)==kInside)
443 {inside = true;}
444 // a very rare case, if a particle remains
445 // on the boundary and escapes the crystal
446 else if (crystallogic->GetSolid()->
447 Inside(xyz0)==kOutside)
448 {inCrystal = false;}//escape the cycle
449 }
450 }
451 while (inCrystal);
452
453 //the angles in the local reference system of the volume
454 tx0 = fCrystalData->AngleXFromLatticeToBox(tx,z);
455 ty0 = ty;
456
457 //set global time
458 fastStep.ProposePrimaryTrackFinalTime(tGlobal);
459 //set final position
461 //set final kinetic energy
463 fastTrack.GetPrimaryTrack()->
464 GetParticleDefinition()->GetPDGMass());
465 //set final momentum direction
466 G4double momentumDirectionZ =
467 1./std::sqrt(1.+std::pow(std::tan(tx0),2)+std::pow(std::tan(ty0),2));
469 G4ThreeVector(momentumDirectionZ*std::tan(tx0),
470 momentumDirectionZ*std::tan(ty0),
471 momentumDirectionZ));
472}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
double z() const
double x() const
double y() const
void ResetRadIntegral()
void GeneratePhoton(G4FastStep &fastStep)
G4double GetParticleNewAngleX()
G4double GetParticleNewTotalEnergy()
G4bool DoRadiation(G4double etotal, G4double mass, G4double angleX, G4double angleY, G4double angleScatteringX, G4double angleScatteringY, G4double step, G4double globalTime, G4ThreeVector coordinateXYZ, G4bool flagEndTrajectory=false)
G4double GetParticleNewAngleY()
G4double GetNewGlobalTime()
const G4ThreeVector & GetParticleNewCoordinateXYZ()
G4ThreeVector CoordinatesFromBoxToLattice(const G4ThreeVector &pos0)
G4double AngleXFromBoxToLattice(G4double tx, G4double z)
G4ThreeVector ChannelChange(G4double &x, G4double &y, G4double &z)
change the channel if necessary, recalculate x o y
G4ThreeVector CoordinatesFromLatticeToBox(const G4ThreeVector &pos)
G4double AngleXFromLatticeToBox(G4double tx, G4double z)
G4double GetLowKineticEnergyLimit(const G4String &particleName)
get cuts
G4double GetLindhardAngleNumberHighLimit(const G4String &particleName)
void SetNumberOfSecondaryTracks(G4int)
void ProposePrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
Definition G4FastStep.cc:96
void ProposePrimaryTrackFinalTime(G4double)
void ProposePrimaryTrackFinalKineticEnergy(G4double)
void ProposePrimaryTrackFinalMomentumDirection(const G4ThreeVector &, G4bool localCoordinates=true)
G4ThreeVector GetPrimaryTrackLocalPosition() const
const G4Track * GetPrimaryTrack() const
G4ThreeVector GetPrimaryTrackLocalDirection() const
G4LogicalVolume * GetEnvelopeLogicalVolume() const
G4VSolid * GetSolid() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetGlobalTime() const
G4double GetTotalEnergy() const
G4double IonizationLosses(G4double dz, G4int ielement)
ionization losses
G4double GetLindhardAngle(G4double etotal, G4double mass)
Calculate the value of the Lindhard angle (!!! the value for a straight crystal)
G4double GetBeta()
get particle velocity/c
void SetGeometryParameters(const G4LogicalVolume *crystallogic)
set geometry parameters from current logical volume
G4double ElectronDensity(G4double x, G4double y)
electron density function
G4double GetSimulationStep(G4double tx, G4double ty)
G4double MinIonizationEnergy(G4double x, G4double y)
minimum energy of ionization function
G4double Ex(G4double x, G4double y)
electric fields produced by crystal lattice
G4ThreeVector CoulombElectronScattering(G4double eMinIonization, G4double electronDensity, G4double step)
multiple and single scattering on electrons
void SetParticleProperties(G4double etotal, G4double mp, G4double charge, G4bool ifhadron)
G4double NuclearDensity(G4double x, G4double y, G4int ielement)
nuclear density function (normalized to average nuclear density)
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68

◆ GetCrystalData()

G4ChannelingFastSimCrystalData * G4ChannelingFastSimModel::GetCrystalData ( )
inline

Definition at line 70 of file G4ChannelingFastSimModel.hh.

70{return fCrystalData;}

◆ GetIfRadiationModelActive()

G4bool G4ChannelingFastSimModel::GetIfRadiationModelActive ( )
inline

Definition at line 74 of file G4ChannelingFastSimModel.hh.

74{return fRad;}

◆ GetLindhardAngleNumberHighLimit() [1/2]

G4double G4ChannelingFastSimModel::GetLindhardAngleNumberHighLimit ( const G4String & particleName)
inline

Definition at line 100 of file G4ChannelingFastSimModel.hh.

101 {return GetLindhardAngleNumberHighLimit(particleTable->
102 FindParticle(particleName)->
103 GetParticleDefinitionID());}

Referenced by DoIt(), GetLindhardAngleNumberHighLimit(), and ModelTrigger().

◆ GetLindhardAngleNumberHighLimit() [2/2]

G4double G4ChannelingFastSimModel::GetLindhardAngleNumberHighLimit ( G4int particleDefinitionID)
inline

Definition at line 109 of file G4ChannelingFastSimModel.hh.

110 {return (fLindhardAngleNumberHighLimit.count(particleDefinitionID) == 1)
111 ? fLindhardAngleNumberHighLimit[particleDefinitionID]
112 : fDefaultLindhardAngleNumberHighLimit;}

◆ GetLowKineticEnergyLimit() [1/2]

G4double G4ChannelingFastSimModel::GetLowKineticEnergyLimit ( const G4String & particleName)
inline

get cuts

Definition at line 96 of file G4ChannelingFastSimModel.hh.

97 {return GetLowKineticEnergyLimit(particleTable->
98 FindParticle(particleName)->
99 GetParticleDefinitionID());}

Referenced by DoIt(), GetLowKineticEnergyLimit(), and ModelTrigger().

◆ GetLowKineticEnergyLimit() [2/2]

G4double G4ChannelingFastSimModel::GetLowKineticEnergyLimit ( G4int particleDefinitionID)
inline

Definition at line 105 of file G4ChannelingFastSimModel.hh.

106 {return (fLowEnergyLimit.count(particleDefinitionID) == 1)
107 ? fLowEnergyLimit[particleDefinitionID]
108 : fDefaultLowEnergyLimit;}

◆ GetMaxPhotonsProducedPerStep()

G4int G4ChannelingFastSimModel::GetMaxPhotonsProducedPerStep ( )
inline

get the maximal number of photons that can be produced per fastStep

Definition at line 115 of file G4ChannelingFastSimModel.hh.

115{return fMaxPhotonsProducedPerStep;}

◆ GetRadiationModel()

G4BaierKatkov * G4ChannelingFastSimModel::GetRadiationModel ( )
inline

Definition at line 72 of file G4ChannelingFastSimModel.hh.

72{return fBaierKatkov;}

◆ Input()

void G4ChannelingFastSimModel::Input ( const G4Material * crystal,
const G4String & lattice )

special functions

Definition at line 476 of file G4ChannelingFastSimModel.cc.

477{
478 //initializing the class with containing all
479 //the crystal material and crystal lattice data and
480 //Channeling scattering and ionization processes
481 fCrystalData = new G4ChannelingFastSimCrystalData();
482 //setting all the crystal material and lattice data
483 fCrystalData->SetMaterialProperties(crystal,lattice);
484
485 //setting default low energy cuts for kinetic energy
486 SetLowKineticEnergyLimit(1*GeV,"proton");
487 SetLowKineticEnergyLimit(1*GeV,"anti_proton");
488 SetLowKineticEnergyLimit(200*MeV,"e-");
489 SetLowKineticEnergyLimit(200*MeV,"e+");
490
491 //set the model high limit of the angle expressed in [Lindhard angle] units
492 SetLindhardAngleNumberHighLimit(100.,"proton");
493 SetLindhardAngleNumberHighLimit(100.,"anti_proton");
496}
void SetMaterialProperties(const G4Material *crystal, const G4String &lattice)
void SetLindhardAngleNumberHighLimit(G4double angleNumber, const G4String &particleName)
void SetLowKineticEnergyLimit(G4double ekinetic, const G4String &particleName)
set cuts

◆ IsApplicable()

G4bool G4ChannelingFastSimModel::IsApplicable ( const G4ParticleDefinition & particleType)
overridevirtual

– IsApplicable

Implements G4VFastSimulationModel.

Definition at line 62 of file G4ChannelingFastSimModel.cc.

63{
64 return std::abs(particleType.GetPDGCharge())>DBL_EPSILON;
65}
#define DBL_EPSILON
Definition templates.hh:66

◆ ModelTrigger()

G4bool G4ChannelingFastSimModel::ModelTrigger ( const G4FastTrack & fastTrack)
overridevirtual

– ModelTrigger

Implements G4VFastSimulationModel.

Definition at line 69 of file G4ChannelingFastSimModel.cc.

70{
71 //default output
72 G4bool modelTrigger = false;
73
74 G4int particleDefinitionID =
76 //kinetic energy
77 G4double ekinetic = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
78
79 //energy cut, at the beginning, to not check everything else
80 if(ekinetic > GetLowKineticEnergyLimit(particleDefinitionID))
81 {
82 //current logical volume
83 G4LogicalVolume* crystallogic = fastTrack.GetEnvelopeLogicalVolume();
84 fCrystalData->SetGeometryParameters(crystallogic);
85
86 G4ThreeVector momentumDirection = fastTrack.GetPrimaryTrackLocalDirection();
87 // the particle angle vs crystal plane or axis
88 G4double angle = std::atan(momentumDirection.x()/momentumDirection.z());
89 //recalculate angle into the lattice reference system
90 angle = fCrystalData->
91 AngleXFromBoxToLattice(angle,
92 (fCrystalData->CoordinatesFromBoxToLattice(
93 fastTrack.GetPrimaryTrackLocalPosition())).z());
94 if (fCrystalData->GetModel()==2)
95 {
96 angle = std::sqrt(angle*angle+
97 std::pow(std::atan(momentumDirection.y()/
98 momentumDirection.z()),2));
99 }
100
101 //particle mass
103 //particle total energy
104 G4double etotal = fastTrack.GetPrimaryTrack()->GetTotalEnergy();
105
106 //Particle position
108 //Step estimate
109 G4double dz0 = fCrystalData->GetMaxSimulationStep(etotal,mass);
110 xyz0 += 2*dz0*momentumDirection;//overestimated particle shift on the next step
111 //in channeling
112
113 //Applies the parameterisation not at the last step, only forward local direction
114 //above low energy limit and below angular limit
115
116 modelTrigger = (crystallogic->GetSolid()->
117 Inside(xyz0)==kInside) &&
118 momentumDirection.z()>0. &&
119 std::abs(angle) <
120 GetLindhardAngleNumberHighLimit(particleDefinitionID) *
121 fCrystalData->GetLindhardAngle(etotal,mass);
122 }
123
124 return modelTrigger;
125}
G4int GetParticleDefinitionID() const
G4double GetKineticEnergy() const
G4double GetMaxSimulationStep(G4double etotal, G4double mass)
Calculate maximal simulation step (standard value for channeling particles)

◆ RadiationModelActivate()

void G4ChannelingFastSimModel::RadiationModelActivate ( )

Definition at line 500 of file G4ChannelingFastSimModel.cc.

501{
502 fRad = true;
503 //activate the Baier-Katkov radiation model
504 fBaierKatkov = new G4BaierKatkov();
505}

◆ SetDefaultLindhardAngleNumberHighLimit()

void G4ChannelingFastSimModel::SetDefaultLindhardAngleNumberHighLimit ( G4double angleNumber)
inline

Definition at line 86 of file G4ChannelingFastSimModel.hh.

87 {fDefaultLindhardAngleNumberHighLimit=angleNumber;}

◆ SetDefaultLowKineticEnergyLimit()

void G4ChannelingFastSimModel::SetDefaultLowKineticEnergyLimit ( G4double ekinetic)
inline

Definition at line 84 of file G4ChannelingFastSimModel.hh.

85 {fDefaultLowEnergyLimit=ekinetic;}

◆ SetLindhardAngleNumberHighLimit()

void G4ChannelingFastSimModel::SetLindhardAngleNumberHighLimit ( G4double angleNumber,
const G4String & particleName )
inline

Definition at line 80 of file G4ChannelingFastSimModel.hh.

81 {fLindhardAngleNumberHighLimit[particleTable->FindParticle(particleName)->
82 GetParticleDefinitionID()]=angleNumber;}
G4ParticleDefinition * FindParticle(G4int PDGEncoding)

Referenced by Input().

◆ SetLowKineticEnergyLimit()

void G4ChannelingFastSimModel::SetLowKineticEnergyLimit ( G4double ekinetic,
const G4String & particleName )
inline

set cuts

Definition at line 77 of file G4ChannelingFastSimModel.hh.

78 {fLowEnergyLimit[particleTable->FindParticle(particleName)->
79 GetParticleDefinitionID()] = ekinetic;}

Referenced by Input().

◆ SetMaxPhotonsProducedPerStep()

void G4ChannelingFastSimModel::SetMaxPhotonsProducedPerStep ( G4double nPhotons)
inline

get the maximal number of photons that can be produced per fastStep Caution: is redundant, if the radiation model is not activated

Definition at line 92 of file G4ChannelingFastSimModel.hh.

93 {fMaxPhotonsProducedPerStep=nPhotons;}

The documentation for this class was generated from the following files: