Geant4 11.4.0
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 ()=default
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 Input (const G4Material *crystal, const G4String &lattice, const G4String &filePath)
void RadiationModelActivate ()
G4ChannelingFastSimCrystalDataGetCrystalData ()
G4BaierKatkovGetRadiationModel ()
G4bool GetIfRadiationModelActive ()
void SetLowKineticEnergyLimit (G4double ekinetic, const G4String &particleName)
 set cuts
void SetLindhardAngleNumberHighLimit (G4double angleNumber, const G4String &particleName)
void SetHighAngleLimit (G4double anglemax, const G4String &particleName)
void SetDefaultLowKineticEnergyLimit (G4double ekinetic)
void SetDefaultLindhardAngleNumberHighLimit (G4double angleNumber)
void SetDefaultHighAngleLimit (G4double anglemax)
void SetMaxPhotonsProducedPerStep (G4double nPhotons)
G4double GetLowKineticEnergyLimit (G4int particleDefinitionID)
 get cuts
G4double GetLindhardAngleNumberHighLimit (G4int particleDefinitionID)
G4double GetHighAngleLimit (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 51 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 ( )
default

Member Function Documentation

◆ DoIt()

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

– User method DoIt

Implements G4VFastSimulationModel.

Definition at line 130 of file G4ChannelingFastSimModel.cc.

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

◆ GetCrystalData()

G4ChannelingFastSimCrystalData * G4ChannelingFastSimModel::GetCrystalData ( )
inline

Definition at line 77 of file G4ChannelingFastSimModel.hh.

77{return fCrystalData;}

Referenced by G4CoherentPairProductionPhysics::ConstructProcess().

◆ GetHighAngleLimit()

G4double G4ChannelingFastSimModel::GetHighAngleLimit ( G4int particleDefinitionID)
inline

Definition at line 116 of file G4ChannelingFastSimModel.hh.

117 {return (fHighAngleLimit.count(particleDefinitionID) == 1)
118 ? fHighAngleLimit[particleDefinitionID]
119 : fDefaultHighAngleLimit;}

Referenced by DoIt(), and ModelTrigger().

◆ GetIfRadiationModelActive()

G4bool G4ChannelingFastSimModel::GetIfRadiationModelActive ( )
inline

Definition at line 81 of file G4ChannelingFastSimModel.hh.

81{return fRad;}

◆ GetLindhardAngleNumberHighLimit()

G4double G4ChannelingFastSimModel::GetLindhardAngleNumberHighLimit ( G4int particleDefinitionID)
inline

Definition at line 112 of file G4ChannelingFastSimModel.hh.

113 {return (fLindhardAngleNumberHighLimit.count(particleDefinitionID) == 1)
114 ? fLindhardAngleNumberHighLimit[particleDefinitionID]
115 : fDefaultLindhardAngleNumberHighLimit;}

Referenced by DoIt(), and ModelTrigger().

◆ GetLowKineticEnergyLimit()

G4double G4ChannelingFastSimModel::GetLowKineticEnergyLimit ( G4int particleDefinitionID)
inline

get cuts

Definition at line 108 of file G4ChannelingFastSimModel.hh.

109 {return (fLowEnergyLimit.count(particleDefinitionID) == 1)
110 ? fLowEnergyLimit[particleDefinitionID]
111 : fDefaultLowEnergyLimit;}

Referenced by DoIt(), and ModelTrigger().

◆ GetMaxPhotonsProducedPerStep()

G4int G4ChannelingFastSimModel::GetMaxPhotonsProducedPerStep ( )
inline

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

Definition at line 122 of file G4ChannelingFastSimModel.hh.

122{return fMaxPhotonsProducedPerStep;}

◆ GetRadiationModel()

G4BaierKatkov * G4ChannelingFastSimModel::GetRadiationModel ( )
inline

Definition at line 79 of file G4ChannelingFastSimModel.hh.

79{return fBaierKatkov;}

◆ Input() [1/2]

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

special functions

Definition at line 67 of file G4ChannelingFastSimModel.hh.

69 {Input(crystal,lattice,"");}
void Input(const G4Material *crystal, const G4String &lattice)
special functions

Referenced by Input().

◆ Input() [2/2]

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

Definition at line 520 of file G4ChannelingFastSimModel.cc.

523{
524 //initializing the class with containing all
525 //the crystal material and crystal lattice data and
526 //Channeling scattering and ionization processes
527 fCrystalData = new G4ChannelingFastSimCrystalData();
528 //setting all the crystal material and lattice data
529 fCrystalData->SetMaterialProperties(crystal,lattice,filePath);
530}

◆ IsApplicable()

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

– IsApplicable

Implements G4VFastSimulationModel.

Definition at line 56 of file G4ChannelingFastSimModel.cc.

57{
58 return std::abs(particleType.GetPDGCharge())>DBL_EPSILON;
59}
#define DBL_EPSILON
Definition templates.hh:66

◆ ModelTrigger()

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

– ModelTrigger

Implements G4VFastSimulationModel.

Definition at line 63 of file G4ChannelingFastSimModel.cc.

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

◆ RadiationModelActivate()

void G4ChannelingFastSimModel::RadiationModelActivate ( )

Definition at line 534 of file G4ChannelingFastSimModel.cc.

535{
536 fRad = true;
537 //activate the Baier-Katkov radiation model
538 fBaierKatkov = new G4BaierKatkov();
539}

◆ SetDefaultHighAngleLimit()

void G4ChannelingFastSimModel::SetDefaultHighAngleLimit ( G4double anglemax)
inline

Definition at line 98 of file G4ChannelingFastSimModel.hh.

99 {fDefaultHighAngleLimit=anglemax;}

◆ SetDefaultLindhardAngleNumberHighLimit()

void G4ChannelingFastSimModel::SetDefaultLindhardAngleNumberHighLimit ( G4double angleNumber)
inline

Definition at line 96 of file G4ChannelingFastSimModel.hh.

97 {fDefaultLindhardAngleNumberHighLimit=angleNumber;}

◆ SetDefaultLowKineticEnergyLimit()

void G4ChannelingFastSimModel::SetDefaultLowKineticEnergyLimit ( G4double ekinetic)
inline

Definition at line 94 of file G4ChannelingFastSimModel.hh.

95 {fDefaultLowEnergyLimit=ekinetic;}

◆ SetHighAngleLimit()

void G4ChannelingFastSimModel::SetHighAngleLimit ( G4double anglemax,
const G4String & particleName )
inline

Definition at line 90 of file G4ChannelingFastSimModel.hh.

91 {fHighAngleLimit[particleTable->FindParticle(particleName)->
92 GetParticleDefinitionID()] = anglemax;}

◆ SetLindhardAngleNumberHighLimit()

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

Definition at line 87 of file G4ChannelingFastSimModel.hh.

88 {fLindhardAngleNumberHighLimit[particleTable->FindParticle(particleName)->
89 GetParticleDefinitionID()]=angleNumber;}

◆ SetLowKineticEnergyLimit()

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

set cuts

Definition at line 84 of file G4ChannelingFastSimModel.hh.

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

◆ 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 104 of file G4ChannelingFastSimModel.hh.

105 {fMaxPhotonsProducedPerStep=nPhotons;}

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