Geant4  10.03
ExExChProcessChanneling.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 //
28 //
29 //
30 
31 #ifndef ExExChProcessChanneling_h
32 #define ExExChProcessChanneling_h 1
33 
34 #include "G4ios.hh"
35 #include "globals.hh"
36 #include "G4VDiscreteProcess.hh"
37 
38 #include "XLatticeManager3.hh"
41 
42 #include "G4VPhysicalVolume.hh"
43 #include "XPhysicalLattice.hh"
45 
47 {
48 public:
49 
50  ExExChProcessChanneling(const G4String& processName =
51  "ExExChProcessChanneling" );
52 
53  virtual ~ExExChProcessChanneling();
54 
55  virtual G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
56 
58 
59  virtual void BuildPhysicsTable(const G4ParticleDefinition&);
60 
61 protected:
62 
63  virtual G4double GetMeanFreePath(const G4Track&,
64  G4double,
66 
67 private:
68 
70 
71 public:
74 
77 
80 
83 
86 
89  };
90 
92  fTransverseVariationMax = aDouble;
93  };
94 
96  void SetTimeStepMin(G4double aDouble) {fTimeStepMin = aDouble;};
97 
99 
100  void SetFileCharacteristicsName(const G4String& vFilename)
101  {fFileCharacteristicsName = vFilename;};
103 
104 private:
105  void UpdateParameters(const G4Track&);
107  void ResetDensity(const G4Track&);
108 
113 
117 
124 
131 
132  G4bool HasLattice(const G4Track&);
136 
141 
143 
144 private:
145  //binding methods
150 
151 private:
152  // hide assignment operator as private
156 
157 private:
159 
163 
166 
168 
169 private:
175 
179 
182 };
183 
184 #endif
185 
186 
187 
188 
189 
190 
191 
192 
193 
194 
void SetFileCharacteristicsName(const G4String &vFilename)
G4ThreeVector ComputeCentrifugalEnergy(const G4Track &, G4ThreeVector)
G4double ComputeCriticalEnergyMinimum(const G4Track &)
XVCrystalCharacteristic * fElectronDensity
G4bool HasLattice(const G4Track &)
G4double ComputePotentialEnergyBent(const G4Track &)
G4double ComputeCriticalRadius(const G4Track &)
CLHEP::Hep3Vector G4ThreeVector
XVCrystalCharacteristic * GetNucleiDensity()
Definition of the XLatticeManager3 class.
G4bool ParticleIsNegative(const G4Track &)
G4ThreeVector ComputeMomentum(const G4Track &, G4StepPoint *)
ExExChProcessChanneling(const G4String &processName="ExExChProcessChanneling")
G4ThreeVector ComputeKineticEnergy(const G4Track &)
Definition of the XVCrystalCharacteristic class.
G4ThreeVector ComputePotentialEnergy(const G4Track &)
G4ThreeVector ComputeTransverseEnergyBent(const G4Track &)
G4ThreeVector ComputeCentrifugalEnergyMaximumVariation(const G4Track &)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
XCrystalIntegratedDensityHub * fIntegratedDensity
Definition of the ExExChParticleUserInfo class.
G4double ComputeCriticalEnergyMaximum(const G4Track &)
G4bool ParticleIsNotOnBoundary(const G4Track &)
XVCrystalCharacteristic * GetElectronDensity()
G4StepPoint * CheckStepPointLatticeForVolume(G4StepPoint *, const G4Track &)
bool G4bool
Definition: G4Types.hh:79
G4bool HasLatticeOnBoundary(const G4Track &)
XVCrystalCharacteristic * GetElectricField()
G4VPhysicalVolume * GetVolume(const G4Track &)
XVCrystalCharacteristic * fNucleiDensity
void SetPotential(XVCrystalCharacteristic *)
void SetTimeStepMin(G4double aDouble)
Definition: G4Step.hh:76
G4double ComputeCriticalEnergyBent(const G4Track &)
G4StepPoint * CheckStepPointLatticeForPosition(G4StepPoint *, const G4Track &)
void ReadFromFileCharacteristics(G4bool)
XVCrystalCharacteristic * fElectricField
void SetIntegratedDensity(XCrystalIntegratedDensityHub *)
G4bool UpdateInitialParameters(const G4Track &)
G4ThreeVector ComputeTransverseEnergy(const G4Track &)
G4ParticleDefinition * GetParticleDefinition(const G4Track &aTrack)
G4bool ParticleIsNotOnBoundaryPost(const G4Track &)
G4ThreeVector ComputePotentialWellCentre(const G4Track &)
void SetElectricField(XVCrystalCharacteristic *)
void UpdateParameters(const G4Track &)
G4double ComputeCriticalEnergyMinimumBent(const G4Track &)
G4bool HasLatticeOnBoundaryPre(const G4Track &)
XVCrystalCharacteristic * GetPotential()
G4double ComputeOscillationPeriod(const G4Track &)
G4double GetChannelingMeanFreePath(const G4Track &)
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
XCrystalIntegratedDensityHub * GetIntegratedDensity()
void SetNucleiDensity(XVCrystalCharacteristic *)
G4bool HasLatticeOnBoundaryPost(const G4Track &)
G4bool UpdateIntegrationStep(const G4Track &, G4ThreeVector &)
ExExChParticleUserInfo * GetInfo(const G4Track &)
void SetTransverseVariationMax(G4double aDouble)
void SetElectronDensity(XVCrystalCharacteristic *)
XPhysicalLattice * GetXPL(const G4Track &)
G4ThreeVector ComputePositionInTheCrystal(G4StepPoint *, const G4Track &)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double ComputeCriticalAngle(const G4Track &)
void ComputeCrystalCharacteristic(const G4Track &)
Definition of the XPhysicalLattice class.
double G4double
Definition: G4Types.hh:76
ExExChProcessChanneling & operator=(const ExExChProcessChanneling &right)
G4ForceCondition
void ResetDensity(const G4Track &)
G4bool ParticleIsNotOnBoundaryPre(const G4Track &)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition of the XCrystalIntegratedDensityHub class.
XVCrystalCharacteristic * fPotentialEnergy