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

#include <G4Channeling.hh>

+ Inheritance diagram for G4Channeling:

Public Member Functions

 G4Channeling ()
 
virtual ~G4Channeling ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aPD)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double GetCriticalAngle (const G4Track &aTrack)
 
G4double GetOscillationPeriod (const G4Track &aTrack)
 
void PosToLattice (G4StepPoint *step, G4ThreeVector &)
 
G4double GetTransverseVariationMax ()
 
void SetTransverseVariationMax (G4double aDouble)
 
G4double GetTimeStepMin ()
 
void SetTimeStepMin (G4double aDouble)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 40 of file G4Channeling.hh.

Constructor & Destructor Documentation

◆ G4Channeling()

G4Channeling::G4Channeling ( )

Definition at line 37 of file G4Channeling.cc.

37 :
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{}
CLHEP::Hep3Vector G4ThreeVector
static G4int GetModelID(const G4int modelIndex)

◆ ~G4Channeling()

G4Channeling::~G4Channeling ( )
virtual

Definition at line 49 of file G4Channeling.cc.

49{;}

Member Function Documentation

◆ BuildPhysicsTable()

virtual void G4Channeling::BuildPhysicsTable ( const G4ParticleDefinition )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 51 of file G4Channeling.hh.

51{;};

◆ GetCriticalAngle()

G4double G4Channeling::GetCriticalAngle ( const G4Track aTrack)
inline

Definition at line 84 of file G4Channeling.hh.

84 {
85 return std::sqrt(2.0*GetMatData(aTrack)->GetPot()->GetMaxMin()
86 /GetPre(aTrack)->GetTotalEnergy());}

Referenced by GetOscillationPeriod().

◆ GetMeanFreePath()

G4double G4Channeling::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 296 of file G4Channeling.cc.

299 {
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}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced
double G4double
Definition: G4Types.hh:83
G4double GetOscillationPeriod(const G4Track &aTrack)
Definition: G4Channeling.hh:87
static G4bool IsLattice(G4LogicalVolume *aLV)
G4VPhysicalVolume * GetVolume() const
G4VPhysicalVolume * GetNextVolume() const
G4LogicalVolume * GetLogicalVolume() const
#define DBL_MAX
Definition: templates.hh:62

◆ GetOscillationPeriod()

G4double G4Channeling::GetOscillationPeriod ( const G4Track aTrack)
inline

Definition at line 87 of file G4Channeling.hh.

87 {
88 return (CLHEP::pi * GetMatData(aTrack)->GetPot()->GetIntSp(0)
89 / GetCriticalAngle(aTrack));
90 }
G4double GetCriticalAngle(const G4Track &aTrack)
Definition: G4Channeling.hh:84

Referenced by GetMeanFreePath().

◆ GetTimeStepMin()

G4double G4Channeling::GetTimeStepMin ( )
inline

Definition at line 117 of file G4Channeling.hh.

117{return fTimeStepMin;};

◆ GetTransverseVariationMax()

G4double G4Channeling::GetTransverseVariationMax ( )
inline

Definition at line 114 of file G4Channeling.hh.

114{return fTransverseVariationMax;};

◆ IsApplicable()

virtual G4bool G4Channeling::IsApplicable ( const G4ParticleDefinition aPD)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 48 of file G4Channeling.hh.

48 {
49 return(aPD.GetPDGCharge() != 0.);
50 };
G4double GetPDGCharge() const

◆ PosToLattice()

void G4Channeling::PosToLattice ( G4StepPoint step,
G4ThreeVector pos 
)

Definition at line 75 of file G4Channeling.cc.

75 {
76 G4TouchableHistory* theTouchable = (G4TouchableHistory*)(step->GetTouchable());
77
78 pos -= theTouchable->GetTranslation();
79 pos = ((*theTouchable->GetRotation()).inverse())(pos);
80}
const G4VTouchable * GetTouchable() const
const G4RotationMatrix * GetRotation(G4int depth=0) const
const G4ThreeVector & GetTranslation(G4int depth=0) const

Referenced by PostStepDoIt().

◆ PostStepDoIt()

G4VParticleChange * G4Channeling::PostStepDoIt ( const G4Track aTrack,
const G4Step  
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 326 of file G4Channeling.cc.

328 {
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}
bool G4bool
Definition: G4Types.hh:86
double z() const
Hep3Vector unit() const
double x() const
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:24
virtual G4ThreeVector GetBR(G4ThreeVector &v3)
void PosToLattice(G4StepPoint *step, G4ThreeVector &)
Definition: G4Channeling.cc:75
const G4ThreeVector & RotateToSolid(G4ThreeVector &dir) const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
void Initialize(const G4Track &) override
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPostStepPoint() const
const G4Step * GetStep() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331

◆ SetTimeStepMin()

void G4Channeling::SetTimeStepMin ( G4double  aDouble)
inline

Definition at line 118 of file G4Channeling.hh.

118{fTimeStepMin = aDouble;};

◆ SetTransverseVariationMax()

void G4Channeling::SetTransverseVariationMax ( G4double  aDouble)
inline

Definition at line 115 of file G4Channeling.hh.

115{fTransverseVariationMax = aDouble;};

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