Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Cerenkov Class Reference

#include <G4Cerenkov.hh>

Inheritance diagram for G4Cerenkov:
Collaboration diagram for G4Cerenkov:

Public Member Functions

 G4Cerenkov (const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
 
 ~G4Cerenkov ()
 
 G4Cerenkov (const G4Cerenkov &right)
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType) override
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
void SetTrackSecondariesFirst (const G4bool state)
 
G4bool GetTrackSecondariesFirst () const
 
void SetMaxBetaChangePerStep (const G4double d)
 
G4double GetMaxBetaChangePerStep () const
 
void SetMaxNumPhotonsPerStep (const G4int NumPhotons)
 
G4int GetMaxNumPhotonsPerStep () const
 
void SetStackPhotons (const G4bool)
 
G4bool GetStackPhotons () const
 
G4int GetNumPhotons () const
 
G4PhysicsTableGetPhysicsTable () const
 
void DumpPhysicsTable () const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
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 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 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
 
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 Attributes

G4PhysicsTablethePhysicsTable
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 81 of file G4Cerenkov.hh.

Constructor & Destructor Documentation

G4Cerenkov::G4Cerenkov ( const G4String processName = "Cerenkov",
G4ProcessType  type = fElectromagnetic 
)
explicit

Definition at line 96 of file G4Cerenkov.cc.

97  : G4VProcess(processName, type),
98  fTrackSecondariesFirst(false),
99  fMaxBetaChange(0.0),
100  fMaxPhotons(0),
101  fStackingFlag(true),
102  fNumPhotons(0)
103 {
105 
106  thePhysicsTable = nullptr;
107 
108  if (verboseLevel>0) {
109  G4cout << GetProcessName() << " is created " << G4endl;
110  }
111 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:212
G4GLOB_DLL std::ostream G4cout
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4Cerenkov::~G4Cerenkov ( )

Definition at line 121 of file G4Cerenkov.cc.

122 {
123  if (thePhysicsTable != nullptr) {
125  delete thePhysicsTable;
126  }
127 }
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:212
void clearAndDestroy()

Here is the call graph for this function:

G4Cerenkov::G4Cerenkov ( const G4Cerenkov right)
explicit

Member Function Documentation

virtual G4VParticleChange* G4Cerenkov::AlongStepDoIt ( const G4Track ,
const G4Step  
)
inlineoverridevirtual

Implements G4VProcess.

Definition at line 148 of file G4Cerenkov.hh.

149  {return nullptr;};
virtual G4double G4Cerenkov::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
)
inlineoverridevirtual

Implements G4VProcess.

Definition at line 133 of file G4Cerenkov.hh.

138  { return -1.0; };
virtual G4VParticleChange* G4Cerenkov::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlineoverridevirtual

Implements G4VProcess.

Definition at line 145 of file G4Cerenkov.hh.

146  {return nullptr;};
virtual G4double G4Cerenkov::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
inlineoverridevirtual

Implements G4VProcess.

Definition at line 140 of file G4Cerenkov.hh.

142  { return -1.0; };
void G4Cerenkov::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 160 of file G4Cerenkov.cc.

161 {
162  if (!thePhysicsTable) BuildThePhysicsTable();
163 }
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:212
void G4Cerenkov::DumpPhysicsTable ( ) const
inline

Definition at line 269 of file G4Cerenkov.hh.

270 {
271  G4int PhysicsTableSize = thePhysicsTable->entries();
273 
274  for (G4int i = 0 ; i < PhysicsTableSize ; i++ ) {
276  v->DumpValues();
277  }
278 }
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:212
int G4int
Definition: G4Types.hh:78
size_t entries() const

Here is the call graph for this function:

G4double G4Cerenkov::GetMaxBetaChangePerStep ( ) const
inline

Definition at line 239 of file G4Cerenkov.hh.

240 {
241  return fMaxBetaChange;
242 }
G4int G4Cerenkov::GetMaxNumPhotonsPerStep ( ) const
inline

Definition at line 245 of file G4Cerenkov.hh.

246 {
247  return fMaxPhotons;
248 }
G4double G4Cerenkov::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition  
)

Definition at line 485 of file G4Cerenkov.cc.

488 {
489  return 1.;
490 }
G4int G4Cerenkov::GetNumPhotons ( ) const
inline

Definition at line 263 of file G4Cerenkov.hh.

264 {
265  return fNumPhotons;
266 }
G4PhysicsTable * G4Cerenkov::GetPhysicsTable ( ) const
inline

Definition at line 281 of file G4Cerenkov.hh.

282 {
283  return thePhysicsTable;
284 }
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:212
G4bool G4Cerenkov::GetStackPhotons ( ) const
inline

Definition at line 257 of file G4Cerenkov.hh.

258 {
259  return fStackingFlag;
260 }
G4bool G4Cerenkov::GetTrackSecondariesFirst ( ) const
inline

Definition at line 233 of file G4Cerenkov.hh.

234 {
235  return fTrackSecondariesFirst;
236 }
G4bool G4Cerenkov::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 133 of file G4Cerenkov.cc.

134 {
135  G4bool result = false;
136 
137  if (aParticleType.GetPDGCharge() != 0.0 &&
138  aParticleType.GetPDGMass() != 0.0 &&
139  aParticleType.GetParticleName() != "chargedgeantino" &&
140  !aParticleType.IsShortLived() ) { result = true; }
141 
142  return result;
143 }
G4double G4ParticleHPJENDLHEData::G4double result
const G4String & GetParticleName() const
bool G4bool
Definition: G4Types.hh:79
G4double GetPDGMass() const
G4double GetPDGCharge() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4VParticleChange * G4Cerenkov::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
overridevirtual

Implements G4VProcess.

Definition at line 169 of file G4Cerenkov.cc.

178 {
180  // Should we ensure that the material is dispersive?
182 
183  aParticleChange.Initialize(aTrack);
184 
185  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
186  const G4Material* aMaterial = aTrack.GetMaterial();
187 
188  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
189  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
190 
191  G4ThreeVector x0 = pPreStepPoint->GetPosition();
192  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
193  G4double t0 = pPreStepPoint->GetGlobalTime();
194 
195  G4MaterialPropertiesTable* aMaterialPropertiesTable =
196  aMaterial->GetMaterialPropertiesTable();
197  if (!aMaterialPropertiesTable) return pParticleChange;
198 
199  G4MaterialPropertyVector* Rindex =
200  aMaterialPropertiesTable->GetProperty("RINDEX");
201  if (!Rindex) return pParticleChange;
202 
203  // particle charge
204  G4double charge = aParticle->GetDefinition()->GetPDGCharge();
205 
206  // particle beta
207  G4double beta = (pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta())*0.5;
208 
209  fNumPhotons = 0;
210 
211  G4double MeanNumberOfPhotons =
212  GetAverageNumberOfPhotons(charge,beta,aMaterial,Rindex);
213 
214  if (MeanNumberOfPhotons <= 0.0) {
215 
216  // return unchanged particle and no secondaries
217 
219 
220  return pParticleChange;
221 
222  }
223 
224  G4double step_length = aStep.GetStepLength();
225 
226  MeanNumberOfPhotons = MeanNumberOfPhotons * step_length;
227 
228  fNumPhotons = (G4int) G4Poisson(MeanNumberOfPhotons);
229 
230  if ( fNumPhotons <= 0 || !fStackingFlag ) {
231 
232  // return unchanged particle and no secondaries
233 
235 
236  return pParticleChange;
237 
238  }
239 
241 
243 
244  if (fTrackSecondariesFirst) {
245  if (aTrack.GetTrackStatus() == fAlive )
247  }
248 
250 
251  G4double Pmin = Rindex->GetMinLowEdgeEnergy();
252  G4double Pmax = Rindex->GetMaxLowEdgeEnergy();
253  G4double dp = Pmax - Pmin;
254 
255  G4double nMax = Rindex->GetMaxValue();
256 
257  G4double BetaInverse = 1./beta;
258 
259  G4double maxCos = BetaInverse / nMax;
260  G4double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
261 
262  G4double beta1 = pPreStepPoint ->GetBeta();
263  G4double beta2 = pPostStepPoint->GetBeta();
264 
265  G4double MeanNumberOfPhotons1 =
266  GetAverageNumberOfPhotons(charge,beta1,aMaterial,Rindex);
267  G4double MeanNumberOfPhotons2 =
268  GetAverageNumberOfPhotons(charge,beta2,aMaterial,Rindex);
269 
270  for (G4int i = 0; i < fNumPhotons; i++) {
271 
272  // Determine photon energy
273 
274  G4double rand;
275  G4double sampledEnergy, sampledRI;
276  G4double cosTheta, sin2Theta;
277 
278  // sample an energy
279 
280  do {
281  rand = G4UniformRand();
282  sampledEnergy = Pmin + rand * dp;
283  sampledRI = Rindex->Value(sampledEnergy);
284  cosTheta = BetaInverse / sampledRI;
285 
286  sin2Theta = (1.0 - cosTheta)*(1.0 + cosTheta);
287  rand = G4UniformRand();
288 
289  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
290  } while (rand*maxSin2 > sin2Theta);
291 
292  // Generate random position of photon on cone surface
293  // defined by Theta
294 
295  rand = G4UniformRand();
296 
297  G4double phi = twopi*rand;
298  G4double sinPhi = std::sin(phi);
299  G4double cosPhi = std::cos(phi);
300 
301  // calculate x,y, and z components of photon energy
302  // (in coord system with primary particle direction
303  // aligned with the z axis)
304 
305  G4double sinTheta = std::sqrt(sin2Theta);
306  G4double px = sinTheta*cosPhi;
307  G4double py = sinTheta*sinPhi;
308  G4double pz = cosTheta;
309 
310  // Create photon momentum direction vector
311  // The momentum direction is still with respect
312  // to the coordinate system where the primary
313  // particle direction is aligned with the z axis
314 
315  G4ParticleMomentum photonMomentum(px, py, pz);
316 
317  // Rotate momentum direction back to global reference
318  // system
319 
320  photonMomentum.rotateUz(p0);
321 
322  // Determine polarization of new photon
323 
324  G4double sx = cosTheta*cosPhi;
325  G4double sy = cosTheta*sinPhi;
326  G4double sz = -sinTheta;
327 
328  G4ThreeVector photonPolarization(sx, sy, sz);
329 
330  // Rotate back to original coord system
331 
332  photonPolarization.rotateUz(p0);
333 
334  // Generate a new photon:
335 
336  G4DynamicParticle* aCerenkovPhoton =
338 
339  aCerenkovPhoton->SetPolarization(photonPolarization.x(),
340  photonPolarization.y(),
341  photonPolarization.z());
342 
343  aCerenkovPhoton->SetKineticEnergy(sampledEnergy);
344 
345  // Generate new G4Track object:
346 
347  G4double NumberOfPhotons, N;
348 
349  do {
350  rand = G4UniformRand();
351  NumberOfPhotons = MeanNumberOfPhotons1 - rand *
352  (MeanNumberOfPhotons1-MeanNumberOfPhotons2);
353  N = G4UniformRand() *
354  std::max(MeanNumberOfPhotons1,MeanNumberOfPhotons2);
355  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
356  } while (N > NumberOfPhotons);
357 
358  G4double delta = rand * aStep.GetStepLength();
359 
360  G4double deltaTime = delta / (pPreStepPoint->GetVelocity()+
361  rand*(pPostStepPoint->GetVelocity()-
362  pPreStepPoint->GetVelocity())*0.5);
363 
364  G4double aSecondaryTime = t0 + deltaTime;
365 
366  G4ThreeVector aSecondaryPosition = x0 + rand * aStep.GetDeltaPosition();
367 
368  G4Track* aSecondaryTrack =
369  new G4Track(aCerenkovPhoton,aSecondaryTime,aSecondaryPosition);
370 
371  aSecondaryTrack->SetTouchableHandle(
373 
374  aSecondaryTrack->SetParentID(aTrack.GetTrackID());
375 
376  aParticleChange.AddSecondary(aSecondaryTrack);
377  }
378 
379  if (verboseLevel>0) {
380  G4cout <<"\n Exiting from G4Cerenkov::DoIt -- NumberOfSecondaries = "
382  }
383 
384  return pParticleChange;
385 }
const int N
Definition: mixmax.h:43
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4int GetNumberOfSecondaries() const
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetStepLength() const
const G4DynamicParticle * GetDynamicParticle() const
G4TrackStatus GetTrackStatus() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4ParticleDefinition * GetDefinition() const
G4double GetVelocity() const
int G4int
Definition: G4Types.hh:78
static constexpr double twopi
Definition: G4SIunits.hh:76
G4StepPoint * GetPreStepPoint() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetPosition() const
G4int GetTrackID() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4double Value(G4double theEnergy, size_t &lastidx) const
void SetKineticEnergy(G4double aEnergy)
G4Material * GetMaterial() const
static G4OpticalPhoton * OpticalPhoton()
virtual void Initialize(const G4Track &)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
void SetNumberOfSecondaries(G4int totSecondaries)
Hep3Vector unit() const
void SetParentID(const G4int aValue)
G4StepPoint * GetPostStepPoint() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void AddSecondary(G4Track *aSecondary)
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double GetPDGCharge() const
G4ThreeVector GetDeltaPosition() const
G4MaterialPropertyVector * GetProperty(const char *key)
const G4TouchableHandle & GetTouchableHandle() const
G4double GetBeta() const

Here is the call graph for this function:

G4double G4Cerenkov::PostStepGetPhysicalInteractionLength ( const G4Track aTrack,
G4double  ,
G4ForceCondition condition 
)
overridevirtual

Implements G4VProcess.

Definition at line 492 of file G4Cerenkov.cc.

496 {
497  *condition = NotForced;
498  G4double StepLimit = DBL_MAX;
499 
500  const G4Material* aMaterial = aTrack.GetMaterial();
501  G4int materialIndex = aMaterial->GetIndex();
502 
503  // If Physics Vector is not defined no Cerenkov photons
504  // this check avoid string comparison below
505 
506  if(!(*thePhysicsTable)[materialIndex]) { return StepLimit; }
507 
508  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
509  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
510 
511  G4double kineticEnergy = aParticle->GetKineticEnergy();
512  const G4ParticleDefinition* particleType = aParticle->GetDefinition();
513  G4double mass = particleType->GetPDGMass();
514 
515  // particle beta
516  G4double beta = aParticle->GetTotalMomentum() /
517  aParticle->GetTotalEnergy();
518  // particle gamma
519  G4double gamma = aParticle->GetTotalEnergy()/mass;
520 
521  G4MaterialPropertiesTable* aMaterialPropertiesTable =
522  aMaterial->GetMaterialPropertiesTable();
523 
524  G4MaterialPropertyVector* Rindex = NULL;
525 
526  if (aMaterialPropertiesTable)
527  Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
528 
529  G4double nMax;
530  if (Rindex) {
531  nMax = Rindex->GetMaxValue();
532  } else {
533  return StepLimit;
534  }
535 
536  G4double BetaMin = 1./nMax;
537  if ( BetaMin >= 1. ) return StepLimit;
538 
539  G4double GammaMin = 1./std::sqrt(1.-BetaMin*BetaMin);
540 
541  if (gamma < GammaMin ) return StepLimit;
542 
543  G4double kinEmin = mass*(GammaMin-1.);
544 
545  G4double RangeMin = G4LossTableManager::Instance()->GetRange(particleType,
546  kinEmin,
547  couple);
548  G4double Range = G4LossTableManager::Instance()->GetRange(particleType,
549  kineticEnergy,
550  couple);
551 
552  G4double Step = Range - RangeMin;
553 // if (Step < 1.*um ) return StepLimit;
554 
555  if (Step > 0. && Step < StepLimit) StepLimit = Step;
556 
557  // If user has defined an average maximum number of photons to
558  // be generated in a Step, then calculate the Step length for
559  // that number of photons.
560 
561  if (fMaxPhotons > 0) {
562 
563  // particle charge
564  const G4double charge = aParticle->GetDefinition()->GetPDGCharge();
565 
566  G4double MeanNumberOfPhotons =
567  GetAverageNumberOfPhotons(charge,beta,aMaterial,Rindex);
568 
569  Step = 0.;
570  if (MeanNumberOfPhotons > 0.0) Step = fMaxPhotons / MeanNumberOfPhotons;
571 
572  if (Step > 0. && Step < StepLimit) StepLimit = Step;
573  }
574 
575  // If user has defined an maximum allowed change in beta per step
576  if (fMaxBetaChange > 0.) {
577 
578  G4double dedx = G4LossTableManager::Instance()->GetDEDX(particleType,
579  kineticEnergy,
580  couple);
581 
582  G4double deltaGamma = gamma - 1./std::sqrt(1.-beta*beta*
583  (1.-fMaxBetaChange)*
584  (1.-fMaxBetaChange));
585 
586  Step = mass * deltaGamma / dedx;
587 
588  if (Step > 0. && Step < StepLimit) StepLimit = Step;
589 
590  }
591 
593  return StepLimit;
594 }
G4double condition(const G4ErrorSymMatrix &m)
static G4LossTableManager * Instance()
G4double GetKineticEnergy() const
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double GetTotalEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
size_t GetIndex() const
Definition: G4Material.hh:262
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4ParticleDefinition * GetDefinition() const
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:212
int G4int
Definition: G4Types.hh:78
G4double GetTotalMomentum() const
G4Material * GetMaterial() const
G4double GetPDGMass() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4MaterialPropertyVector * GetProperty(const char *key)

Here is the call graph for this function:

void G4Cerenkov::SetMaxBetaChangePerStep ( const G4double  d)

Definition at line 150 of file G4Cerenkov.cc.

151 {
152  fMaxBetaChange = value*CLHEP::perCent;
153 }
const XML_Char int const XML_Char * value
Definition: expat.h:331
static constexpr double perCent

Here is the caller graph for this function:

void G4Cerenkov::SetMaxNumPhotonsPerStep ( const G4int  NumPhotons)

Definition at line 155 of file G4Cerenkov.cc.

156 {
157  fMaxPhotons = NumPhotons;
158 }

Here is the caller graph for this function:

void G4Cerenkov::SetStackPhotons ( const G4bool  stackingFlag)
inline

Definition at line 251 of file G4Cerenkov.hh.

252 {
253  fStackingFlag = stackingFlag;
254 }

Here is the caller graph for this function:

void G4Cerenkov::SetTrackSecondariesFirst ( const G4bool  state)

Definition at line 145 of file G4Cerenkov.cc.

146 {
147  fTrackSecondariesFirst = state;
148 }

Here is the caller graph for this function:

Member Data Documentation

G4PhysicsTable* G4Cerenkov::thePhysicsTable
protected

Definition at line 212 of file G4Cerenkov.hh.


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