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

#include <G4Scintillation.hh>

Inheritance diagram for G4Scintillation:
Collaboration diagram for G4Scintillation:

Public Member Functions

 G4Scintillation (const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic)
 
 ~G4Scintillation ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType) override
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
G4double GetScintillationYieldByParticleType (const G4Track &aTrack, const G4Step &aStep)
 
void SetTrackSecondariesFirst (const G4bool state)
 
G4bool GetTrackSecondariesFirst () const
 
void SetFiniteRiseTime (const G4bool state)
 
G4bool GetFiniteRiseTime () const
 
void SetScintillationYieldFactor (const G4double yieldfactor)
 
G4double GetScintillationYieldFactor () const
 
void SetScintillationExcitationRatio (const G4double ratio)
 
G4double GetScintillationExcitationRatio () const
 
G4PhysicsTableGetFastIntegralTable () const
 
G4PhysicsTableGetSlowIntegralTable () const
 
void AddSaturation (G4EmSaturation *sat)
 
void RemoveSaturation ()
 
G4EmSaturationGetSaturation () const
 
void SetScintillationByParticleType (const G4bool)
 
G4bool GetScintillationByParticleType () const
 
void SetScintillationTrackInfo (const G4bool trackType)
 
G4bool GetScintillationTrackInfo () const
 
void SetStackPhotons (const G4bool)
 
G4bool GetStackPhotons () const
 
G4int GetNumPhotons () const
 
void DumpPhysicsTable () const
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- 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 Member Functions

void BuildThePhysicsTable ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4PhysicsTablefFastIntegralTable
 
G4PhysicsTablefSlowIntegralTable
 
- 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)
 

Detailed Description

Definition at line 88 of file G4Scintillation.hh.

Constructor & Destructor Documentation

G4Scintillation::G4Scintillation ( const G4String processName = "Scintillation",
G4ProcessType  type = fElectromagnetic 
)
explicit

Definition at line 111 of file G4Scintillation.cc.

113  : G4VRestDiscreteProcess(processName, type) ,
114  fTrackSecondariesFirst(false),
115  fFiniteRiseTime(false),
116  fYieldFactor(1.0),
117  fExcitationRatio(1.0),
118  fScintillationByParticleType(false),
119  fScintillationTrackInfo(false),
120  fStackingFlag(true),
121  fNumPhotons(0),
122  fEmSaturation(nullptr)
123 {
125 
126 #ifdef G4DEBUG_SCINTILLATION
127  ScintTrackEDep = 0.;
128  ScintTrackYield = 0.;
129 #endif
130 
131  fFastIntegralTable = nullptr;
132  fSlowIntegralTable = nullptr;
133 
134  if (verboseLevel>0) {
135  G4cout << GetProcessName() << " is created " << G4endl;
136  }
137 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4PhysicsTable * fFastIntegralTable
G4PhysicsTable * fSlowIntegralTable
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:

G4Scintillation::~G4Scintillation ( )

Definition at line 143 of file G4Scintillation.cc.

144 {
145  if (fFastIntegralTable != nullptr) {
147  delete fFastIntegralTable;
148  }
149  if (fSlowIntegralTable != nullptr) {
151  delete fSlowIntegralTable;
152  }
153 }
G4PhysicsTable * fFastIntegralTable
G4PhysicsTable * fSlowIntegralTable
void clearAndDestroy()

Here is the call graph for this function:

Member Function Documentation

void G4Scintillation::AddSaturation ( G4EmSaturation sat)
inline

Definition at line 349 of file G4Scintillation.hh.

350 {
351  fEmSaturation = sat;
352 }

Here is the caller graph for this function:

G4VParticleChange * G4Scintillation::AtRestDoIt ( const G4Track aTrack,
const G4Step aStep 
)
overridevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 179 of file G4Scintillation.cc.

184 {
185  return G4Scintillation::PostStepDoIt(aTrack, aStep);
186 }
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override

Here is the call graph for this function:

void G4Scintillation::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 159 of file G4Scintillation.cc.

160 {
161  if (fFastIntegralTable != nullptr) {
163  delete fFastIntegralTable;
164  fFastIntegralTable = nullptr;
165  }
166  if (fSlowIntegralTable != nullptr) {
168  delete fSlowIntegralTable;
169  fSlowIntegralTable = nullptr;
170  }
172 }
G4PhysicsTable * fFastIntegralTable
G4PhysicsTable * fSlowIntegralTable
void clearAndDestroy()

Here is the call graph for this function:

void G4Scintillation::BuildThePhysicsTable ( )
protected

Definition at line 501 of file G4Scintillation.cc.

502 {
503  if (fFastIntegralTable && fSlowIntegralTable) return;
504 
505  const G4MaterialTable* theMaterialTable =
507  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
508 
509  // create new physics table
510 
512  new G4PhysicsTable(numOfMaterials);
514  new G4PhysicsTable(numOfMaterials);
515 
516  // loop for materials
517 
518  for (G4int i=0 ; i < numOfMaterials; i++)
519  {
520  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
522  G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector =
524 
525  // Retrieve vector of scintillation wavelength intensity for
526  // the material from the material's optical properties table.
527 
528  G4Material* aMaterial = (*theMaterialTable)[i];
529 
530  G4MaterialPropertiesTable* aMaterialPropertiesTable =
531  aMaterial->GetMaterialPropertiesTable();
532 
533  if (aMaterialPropertiesTable) {
534 
535  G4MaterialPropertyVector* theFastLightVector =
536  aMaterialPropertiesTable->GetProperty("FASTCOMPONENT");
537 
538  if (theFastLightVector) {
539 
540  // Retrieve the first intensity point in vector
541  // of (photon energy, intensity) pairs
542 
543  G4double currentIN = (*theFastLightVector)[0];
544 
545  if (currentIN >= 0.0) {
546 
547  // Create first (photon energy, Scintillation
548  // Integral pair
549 
550  G4double currentPM = theFastLightVector->Energy(0);
551 
552  G4double currentCII = 0.0;
553 
554  aPhysicsOrderedFreeVector->
555  InsertValues(currentPM , currentCII);
556 
557  // Set previous values to current ones prior to loop
558 
559  G4double prevPM = currentPM;
560  G4double prevCII = currentCII;
561  G4double prevIN = currentIN;
562 
563  // loop over all (photon energy, intensity)
564  // pairs stored for this material
565 
566  for (size_t ii = 1;
567  ii < theFastLightVector->GetVectorLength();
568  ++ii)
569  {
570  currentPM = theFastLightVector->Energy(ii);
571  currentIN = (*theFastLightVector)[ii];
572 
573  currentCII = 0.5 * (prevIN + currentIN);
574 
575  currentCII = prevCII +
576  (currentPM - prevPM) * currentCII;
577 
578  aPhysicsOrderedFreeVector->
579  InsertValues(currentPM, currentCII);
580 
581  prevPM = currentPM;
582  prevCII = currentCII;
583  prevIN = currentIN;
584  }
585 
586  }
587  }
588 
589  G4MaterialPropertyVector* theSlowLightVector =
590  aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
591 
592  if (theSlowLightVector) {
593 
594  // Retrieve the first intensity point in vector
595  // of (photon energy, intensity) pairs
596 
597  G4double currentIN = (*theSlowLightVector)[0];
598 
599  if (currentIN >= 0.0) {
600 
601  // Create first (photon energy, Scintillation
602  // Integral pair
603 
604  G4double currentPM = theSlowLightVector->Energy(0);
605 
606  G4double currentCII = 0.0;
607 
608  bPhysicsOrderedFreeVector->
609  InsertValues(currentPM , currentCII);
610 
611  // Set previous values to current ones prior to loop
612 
613  G4double prevPM = currentPM;
614  G4double prevCII = currentCII;
615  G4double prevIN = currentIN;
616 
617  // loop over all (photon energy, intensity)
618  // pairs stored for this material
619 
620  for (size_t ii = 1;
621  ii < theSlowLightVector->GetVectorLength();
622  ++ii)
623  {
624  currentPM = theSlowLightVector->Energy(ii);
625  currentIN = (*theSlowLightVector)[ii];
626 
627  currentCII = 0.5 * (prevIN + currentIN);
628 
629  currentCII = prevCII +
630  (currentPM - prevPM) * currentCII;
631 
632  bPhysicsOrderedFreeVector->
633  InsertValues(currentPM, currentCII);
634 
635  prevPM = currentPM;
636  prevCII = currentCII;
637  prevIN = currentIN;
638  }
639 
640  }
641  }
642  }
643 
644  // The scintillation integral(s) for a given material
645  // will be inserted in the table(s) according to the
646  // position of the material in the material table.
647 
648  fFastIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
649  fSlowIntegralTable->insertAt(i,bPhysicsOrderedFreeVector);
650 
651  }
652 }
G4PhysicsTable * fFastIntegralTable
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
G4PhysicsTable * fSlowIntegralTable
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
G4double Energy(size_t index) const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:594
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
void insertAt(size_t, G4PhysicsVector *)
double G4double
Definition: G4Types.hh:76
G4MaterialPropertyVector * GetProperty(const char *key)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scintillation::DumpPhysicsTable ( ) const
inline

Definition at line 403 of file G4Scintillation.hh.

404 {
405  if (fFastIntegralTable) {
406  G4int PhysicsTableSize = fFastIntegralTable->entries();
408 
409  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
410  {
412  v->DumpValues();
413  }
414  }
415 
416  if (fSlowIntegralTable) {
417  G4int PhysicsTableSize = fSlowIntegralTable->entries();
419 
420  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
421  {
423  v->DumpValues();
424  }
425  }
426 }
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const
G4PhysicsTable * fFastIntegralTable
G4PhysicsTable * fSlowIntegralTable
int G4int
Definition: G4Types.hh:78
size_t entries() const

Here is the call graph for this function:

G4PhysicsTable * G4Scintillation::GetFastIntegralTable ( ) const
inline

Definition at line 343 of file G4Scintillation.hh.

344 {
345  return fFastIntegralTable;
346 }
G4PhysicsTable * fFastIntegralTable
G4bool G4Scintillation::GetFiniteRiseTime ( ) const
inline

Definition at line 307 of file G4Scintillation.hh.

308 {
309  return fFiniteRiseTime;
310 }
G4double G4Scintillation::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition condition 
)
overridevirtual

Implements G4VRestDiscreteProcess.

Definition at line 668 of file G4Scintillation.cc.

671 {
673 
674  return DBL_MAX;
675 
676 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
G4double G4Scintillation::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition condition 
)
overridevirtual

Implements G4VRestDiscreteProcess.

Definition at line 682 of file G4Scintillation.cc.

684 {
685  *condition = Forced;
686 
687  return DBL_MAX;
688 
689 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
G4int G4Scintillation::GetNumPhotons ( ) const
inline

Definition at line 397 of file G4Scintillation.hh.

398 {
399  return fNumPhotons;
400 }
G4EmSaturation * G4Scintillation::GetSaturation ( ) const
inline

Definition at line 361 of file G4Scintillation.hh.

362 {
363  return fEmSaturation;
364 }
G4bool G4Scintillation::GetScintillationByParticleType ( ) const
inline

Definition at line 367 of file G4Scintillation.hh.

368 {
369  return fScintillationByParticleType;
370 }
G4double G4Scintillation::GetScintillationExcitationRatio ( ) const
inline

Definition at line 331 of file G4Scintillation.hh.

332 {
333  return fExcitationRatio;
334 }
G4bool G4Scintillation::GetScintillationTrackInfo ( ) const
inline

Definition at line 379 of file G4Scintillation.hh.

380 {
381  return fScintillationTrackInfo;
382 }
G4double G4Scintillation::GetScintillationYieldByParticleType ( const G4Track aTrack,
const G4Step aStep 
)

Definition at line 714 of file G4Scintillation.cc.

715 {
717  // Get the scintillation yield vector //
719 
721 
722  G4MaterialPropertyVector *Scint_Yield_Vector = nullptr;
723 
724  G4MaterialPropertiesTable *aMaterialPropertiesTable
726 
727  // Get the G4MaterialPropertyVector containing the scintillation
728  // yield as a function of the energy deposited and particle type
729 
730  // Protons
731  if(pDef==G4Proton::ProtonDefinition())
732  Scint_Yield_Vector = aMaterialPropertiesTable->
733  GetProperty("PROTONSCINTILLATIONYIELD");
734 
735  // Deuterons
736  else if(pDef==G4Deuteron::DeuteronDefinition())
737  Scint_Yield_Vector = aMaterialPropertiesTable->
738  GetProperty("DEUTERONSCINTILLATIONYIELD");
739 
740  // Tritons
741  else if(pDef==G4Triton::TritonDefinition())
742  Scint_Yield_Vector = aMaterialPropertiesTable->
743  GetProperty("TRITONSCINTILLATIONYIELD");
744 
745  // Alphas
746  else if(pDef==G4Alpha::AlphaDefinition())
747  Scint_Yield_Vector = aMaterialPropertiesTable->
748  GetProperty("ALPHASCINTILLATIONYIELD");
749 
750  // Ions (particles derived from G4VIon and G4Ions) and recoil ions
751  // below the production cut from neutrons after hElastic
752  else if(pDef->GetParticleType()== "nucleus" ||
754  Scint_Yield_Vector = aMaterialPropertiesTable->
755  GetProperty("IONSCINTILLATIONYIELD");
756 
757  // Electrons (must also account for shell-binding energy
758  // attributed to gamma from standard photoelectric effect)
759  else if(pDef==G4Electron::ElectronDefinition() ||
760  pDef==G4Gamma::GammaDefinition())
761  Scint_Yield_Vector = aMaterialPropertiesTable->
762  GetProperty("ELECTRONSCINTILLATIONYIELD");
763 
764  // Default for particles not enumerated/listed above
765  else
766  Scint_Yield_Vector = aMaterialPropertiesTable->
767  GetProperty("ELECTRONSCINTILLATIONYIELD");
768 
769  // If the user has specified none of the above particles then the
770  // default is the electron scintillation yield
771  if(!Scint_Yield_Vector)
772  Scint_Yield_Vector = aMaterialPropertiesTable->
773  GetProperty("ELECTRONSCINTILLATIONYIELD");
774 
775  // Throw an exception if no scintillation yield vector is found
776  if (!Scint_Yield_Vector) {
778  ed << "\nG4Scintillation::PostStepDoIt(): "
779  << "Request for scintillation yield for energy deposit and particle\n"
780  << "type without correct entry in MaterialPropertiesTable.\n"
781  << "ScintillationByParticleType requires at minimum that \n"
782  << "ELECTRONSCINTILLATIONYIELD is set by the user\n"
783  << G4endl;
784  G4String comments = "Missing MaterialPropertiesTable entry - No correct entry in MaterialPropertiesTable";
785  G4Exception("G4Scintillation::PostStepDoIt","Scint01",
786  FatalException,ed,comments);
787  }
788 
790  // Calculate the scintillation light //
792  // To account for potential nonlinearity and scintillation photon
793  // density along the track, light (L) is produced according to:
794  //
795  // L_currentStep = L(PreStepKE) - L(PreStepKE - EDep)
796 
797  G4double ScintillationYield = 0.;
798 
799  G4double StepEnergyDeposit = aStep.GetTotalEnergyDeposit();
800  G4double PreStepKineticEnergy = aStep.GetPreStepPoint()->GetKineticEnergy();
801 
802  if(PreStepKineticEnergy <= Scint_Yield_Vector->GetMaxEnergy()){
803  G4double Yield1 = Scint_Yield_Vector->Value(PreStepKineticEnergy);
804  G4double Yield2 = Scint_Yield_Vector->
805  Value(PreStepKineticEnergy - StepEnergyDeposit);
806  ScintillationYield = Yield1 - Yield2;
807  } else {
809  ed << "\nG4Scintillation::GetScintillationYieldByParticleType(): Request\n"
810  << "for scintillation light yield above the available energy range\n"
811  << "specifed in G4MaterialPropertiesTable. A linear interpolation\n"
812  << "will be performed to compute the scintillation light yield using\n"
813  << "(L_max / E_max) as the photon yield per unit energy."
814  << G4endl;
815  G4String cmt = "\nScintillation yield may be unphysical!\n";
816  G4Exception("G4Scintillation::GetScintillationYieldByParticleType()",
817  "Scint03", JustWarning, ed, cmt);
818 
819  G4double LinearYield = Scint_Yield_Vector->GetMaxValue()
820  / Scint_Yield_Vector->GetMaxEnergy();
821 
822  // Units: [# scintillation photons]
823  ScintillationYield = LinearYield * StepEnergyDeposit;
824  }
825 
826 #ifdef G4DEBUG_SCINTILLATION
827 
828  // Increment track aggregators
829  ScintTrackYield += ScintillationYield;
830  ScintTrackEDep += StepEnergyDeposit;
831 
832  G4cout << "\n--- G4Scintillation::GetScintillationYieldByParticleType() ---\n"
833  << "--\n"
834  << "-- Name = " << aTrack.GetParticleDefinition()->GetParticleName() << "\n"
835  << "-- TrackID = " << aTrack.GetTrackID() << "\n"
836  << "-- ParentID = " << aTrack.GetParentID() << "\n"
837  << "-- Current KE = " << aTrack.GetKineticEnergy()/MeV << " MeV\n"
838  << "-- Step EDep = " << aStep.GetTotalEnergyDeposit()/MeV << " MeV\n"
839  << "-- Track EDep = " << ScintTrackEDep/MeV << " MeV\n"
840  << "-- Vertex KE = " << aTrack.GetVertexKineticEnergy()/MeV << " MeV\n"
841  << "-- Step yield = " << ScintillationYield << " photons\n"
842  << "-- Track yield = " << ScintTrackYield << " photons\n"
843  << G4endl;
844 
845  // The track has terminated within or has left the scintillator volume
846  if( (aTrack.GetTrackStatus() == fStopButAlive) or
847  (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary) ){
848 
849  // Reset aggregators for the next track
850  ScintTrackEDep = 0.;
851  ScintTrackYield = 0.;
852  }
853 
854 #endif
855 
856  return ScintillationYield;
857 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Triton * TritonDefinition()
Definition: G4Triton.cc:90
G4int GetParentID() const
G4double GetMaxEnergy() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4DynamicParticle * GetDynamicParticle() const
G4StepStatus GetStepStatus() const
G4TrackStatus GetTrackStatus() const
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4ParticleDefinition * GetDefinition() const
const G4String & GetParticleName() const
G4StepPoint * GetPreStepPoint() const
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4double GetVertexKineticEnergy() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4String & GetParticleType() const
G4int GetTrackID() const
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double GetTotalEnergyDeposit() const
G4Material * GetMaterial() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
G4StepPoint * GetPostStepPoint() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:89
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:84
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4Scintillation::GetScintillationYieldFactor ( ) const
inline

Definition at line 319 of file G4Scintillation.hh.

320 {
321  return fYieldFactor;
322 }
G4PhysicsTable * G4Scintillation::GetSlowIntegralTable ( ) const
inline

Definition at line 337 of file G4Scintillation.hh.

338 {
339  return fSlowIntegralTable;
340 }
G4PhysicsTable * fSlowIntegralTable
G4bool G4Scintillation::GetStackPhotons ( ) const
inline

Definition at line 391 of file G4Scintillation.hh.

392 {
393  return fStackingFlag;
394 }
G4bool G4Scintillation::GetTrackSecondariesFirst ( ) const
inline

Definition at line 295 of file G4Scintillation.hh.

296 {
297  return fTrackSecondariesFirst;
298 }
G4bool G4Scintillation::IsApplicable ( const G4ParticleDefinition aParticleType)
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 280 of file G4Scintillation.hh.

281 {
282  if (aParticleType.GetParticleName() == "opticalphoton") return false;
283  if (aParticleType.IsShortLived()) return false;
284 
285  return true;
286 }
const G4String & GetParticleName() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 192 of file G4Scintillation.cc.

199 {
200  aParticleChange.Initialize(aTrack);
201 
202  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
203  const G4Material* aMaterial = aTrack.GetMaterial();
204 
205  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
206  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
207 
208  G4ThreeVector x0 = pPreStepPoint->GetPosition();
209  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
210  G4double t0 = pPreStepPoint->GetGlobalTime();
211 
212  G4double TotalEnergyDeposit = aStep.GetTotalEnergyDeposit();
213 
214  G4MaterialPropertiesTable* aMaterialPropertiesTable =
215  aMaterial->GetMaterialPropertiesTable();
216  if (!aMaterialPropertiesTable)
217  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
218 
219  G4MaterialPropertyVector* Fast_Intensity =
220  aMaterialPropertiesTable->GetProperty("FASTCOMPONENT");
221  G4MaterialPropertyVector* Slow_Intensity =
222  aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
223 
224  if (!Fast_Intensity && !Slow_Intensity )
225  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
226 
227  G4int nscnt = 1;
228  if (Fast_Intensity && Slow_Intensity) nscnt = 2;
229 
230  G4double ScintillationYield = 0.;
231 
232  // Scintillation depends on particle type, energy deposited
233  if (fScintillationByParticleType) {
234 
235  ScintillationYield =
237 
238  // The default linear scintillation process
239  } else {
240 
241  ScintillationYield = aMaterialPropertiesTable->
242  GetConstProperty("SCINTILLATIONYIELD");
243 
244  // Units: [# scintillation photons / MeV]
245  ScintillationYield *= fYieldFactor;
246  }
247 
248  G4double ResolutionScale = aMaterialPropertiesTable->
249  GetConstProperty("RESOLUTIONSCALE");
250 
251  // Birks law saturation:
252 
253  //G4double constBirks = 0.0;
254 
255  //constBirks = aMaterial->GetIonisation()->GetBirksConstant();
256 
257  G4double MeanNumberOfPhotons;
258 
259  // Birk's correction via fEmSaturation and specifying scintillation by
260  // by particle type are physically mutually exclusive
261 
262  if (fScintillationByParticleType)
263  MeanNumberOfPhotons = ScintillationYield;
264  else if (fEmSaturation)
265  MeanNumberOfPhotons = ScintillationYield*
266  (fEmSaturation->VisibleEnergyDepositionAtAStep(&aStep));
267  else
268  MeanNumberOfPhotons = ScintillationYield*TotalEnergyDeposit;
269 
270  if (MeanNumberOfPhotons > 10.)
271  {
272  G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
273  fNumPhotons=G4int(G4RandGauss::shoot(MeanNumberOfPhotons,sigma)+0.5);
274  }
275  else
276  {
277  fNumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
278  }
279 
280  if ( fNumPhotons <= 0 || !fStackingFlag )
281  {
282  // return unchanged particle and no secondaries
283 
285 
286  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
287  }
288 
290 
292 
293  if (fTrackSecondariesFirst) {
294  if (aTrack.GetTrackStatus() == fAlive )
296  }
297 
299 
300  G4int materialIndex = aMaterial->GetIndex();
301 
302  // Retrieve the Scintillation Integral for this material
303  // new G4PhysicsOrderedFreeVector allocated to hold CII's
304 
305  G4int Num = fNumPhotons;
306 
307  for (G4int scnt = 1; scnt <= nscnt; scnt++) {
308 
309  G4double ScintillationTime = 0.*ns;
310  G4double ScintillationRiseTime = 0.*ns;
311  G4PhysicsOrderedFreeVector* ScintillationIntegral = nullptr;
312  G4ScintillationType ScintillationType = Slow;
313 
314  if (scnt == 1) {
315  if (nscnt == 1) {
316  if (Fast_Intensity) {
317  ScintillationTime = aMaterialPropertiesTable->
318  GetConstProperty("FASTTIMECONSTANT");
319  if (fFiniteRiseTime) {
320  ScintillationRiseTime = aMaterialPropertiesTable->
321  GetConstProperty("FASTSCINTILLATIONRISETIME");
322  }
323  ScintillationType = Fast;
324  ScintillationIntegral =
326  ((*fFastIntegralTable)(materialIndex));
327  }
328  if (Slow_Intensity) {
329  ScintillationTime = aMaterialPropertiesTable->
330  GetConstProperty("SLOWTIMECONSTANT");
331  if (fFiniteRiseTime) {
332  ScintillationRiseTime = aMaterialPropertiesTable->
333  GetConstProperty("SLOWSCINTILLATIONRISETIME");
334  }
335  ScintillationType = Slow;
336  ScintillationIntegral =
338  ((*fSlowIntegralTable)(materialIndex));
339  }
340  }
341  else {
342  G4double yieldRatio = aMaterialPropertiesTable->
343  GetConstProperty("YIELDRATIO");
344  if ( fExcitationRatio == 1.0 || fExcitationRatio == 0.0) {
345  Num = G4int (std::min(yieldRatio,1.0) * fNumPhotons);
346  }
347  else {
348  Num = G4int (std::min(fExcitationRatio,1.0) * fNumPhotons);
349  }
350  ScintillationTime = aMaterialPropertiesTable->
351  GetConstProperty("FASTTIMECONSTANT");
352  if (fFiniteRiseTime) {
353  ScintillationRiseTime = aMaterialPropertiesTable->
354  GetConstProperty("FASTSCINTILLATIONRISETIME");
355  }
356  ScintillationType = Fast;
357  ScintillationIntegral =
359  ((*fFastIntegralTable)(materialIndex));
360  }
361  }
362  else {
363  Num = fNumPhotons - Num;
364  ScintillationTime = aMaterialPropertiesTable->
365  GetConstProperty("SLOWTIMECONSTANT");
366  if (fFiniteRiseTime) {
367  ScintillationRiseTime = aMaterialPropertiesTable->
368  GetConstProperty("SLOWSCINTILLATIONRISETIME");
369  }
370  ScintillationType = Slow;
371  ScintillationIntegral =
373  ((*fSlowIntegralTable)(materialIndex));
374  }
375 
376  if (!ScintillationIntegral) continue;
377 
378  // Max Scintillation Integral
379 
380  G4double CIImax = ScintillationIntegral->GetMaxValue();
381 
382  for (G4int i = 0; i < Num; i++) {
383 
384  // Determine photon energy
385 
386  G4double CIIvalue = G4UniformRand()*CIImax;
387  G4double sampledEnergy =
388  ScintillationIntegral->GetEnergy(CIIvalue);
389 
390  if (verboseLevel>1) {
391  G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
392  G4cout << "CIIvalue = " << CIIvalue << G4endl;
393  }
394 
395  // Generate random photon direction
396 
397  G4double cost = 1. - 2.*G4UniformRand();
398  G4double sint = std::sqrt((1.-cost)*(1.+cost));
399 
400  G4double phi = twopi*G4UniformRand();
401  G4double sinp = std::sin(phi);
402  G4double cosp = std::cos(phi);
403 
404  G4double px = sint*cosp;
405  G4double py = sint*sinp;
406  G4double pz = cost;
407 
408  // Create photon momentum direction vector
409 
410  G4ParticleMomentum photonMomentum(px, py, pz);
411 
412  // Determine polarization of new photon
413 
414  G4double sx = cost*cosp;
415  G4double sy = cost*sinp;
416  G4double sz = -sint;
417 
418  G4ThreeVector photonPolarization(sx, sy, sz);
419 
420  G4ThreeVector perp = photonMomentum.cross(photonPolarization);
421 
422  phi = twopi*G4UniformRand();
423  sinp = std::sin(phi);
424  cosp = std::cos(phi);
425 
426  photonPolarization = cosp * photonPolarization + sinp * perp;
427 
428  photonPolarization = photonPolarization.unit();
429 
430  // Generate a new photon:
431 
432  G4DynamicParticle* aScintillationPhoton =
434  photonMomentum);
435  aScintillationPhoton->SetPolarization
436  (photonPolarization.x(),
437  photonPolarization.y(),
438  photonPolarization.z());
439 
440  aScintillationPhoton->SetKineticEnergy(sampledEnergy);
441 
442  // Generate new G4Track object:
443 
444  G4double rand;
445 
446  if (aParticle->GetDefinition()->GetPDGCharge() != 0) {
447  rand = G4UniformRand();
448  } else {
449  rand = 1.0;
450  }
451 
452  G4double delta = rand * aStep.GetStepLength();
453  G4double deltaTime = delta / (pPreStepPoint->GetVelocity()+
454  rand*(pPostStepPoint->GetVelocity()-
455  pPreStepPoint->GetVelocity())/2.);
456 
457  // emission time distribution
458  if (ScintillationRiseTime==0.0) {
459  deltaTime = deltaTime -
460  ScintillationTime * std::log( G4UniformRand() );
461  } else {
462  deltaTime = deltaTime +
463  sample_time(ScintillationRiseTime, ScintillationTime);
464  }
465 
466  G4double aSecondaryTime = t0 + deltaTime;
467 
468  G4ThreeVector aSecondaryPosition =
469  x0 + rand * aStep.GetDeltaPosition();
470 
471  G4Track* aSecondaryTrack = new G4Track(aScintillationPhoton,
472  aSecondaryTime,
473  aSecondaryPosition);
474 
475  aSecondaryTrack->SetTouchableHandle(
477  // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
478 
479  aSecondaryTrack->SetParentID(aTrack.GetTrackID());
480 
481  if (fScintillationTrackInfo) aSecondaryTrack->
482  SetUserInformation(new G4ScintillationTrackInformation(ScintillationType));
483 
484  aParticleChange.AddSecondary(aSecondaryTrack);
485 
486  }
487  }
488 
489  if (verboseLevel>0) {
490  G4cout << "\n Exiting from G4Scintillation::DoIt -- NumberOfSecondaries = "
492  }
493 
494  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
495 }
ThreeVector shoot(const G4int Ap, const G4int Af)
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4int GetNumberOfSecondaries() const
G4int verboseLevel
Definition: G4VProcess.hh:368
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double GetStepLength() const
const G4DynamicParticle * GetDynamicParticle() const
size_t GetIndex() const
Definition: G4Material.hh:262
G4PhysicsTable * fFastIntegralTable
G4TrackStatus GetTrackStatus() const
G4PhysicsTable * fSlowIntegralTable
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 VisibleEnergyDepositionAtAStep(const G4Step *) const
void SetKineticEnergy(G4double aEnergy)
G4double GetEnergy(G4double aValue)
G4double GetTotalEnergyDeposit() const
G4Material * GetMaterial() const
static G4OpticalPhoton * OpticalPhoton()
virtual void Initialize(const G4Track &)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
void SetNumberOfSecondaries(G4int totSecondaries)
Hep3Vector unit() const
void SetParentID(const G4int aValue)
G4StepPoint * GetPostStepPoint() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void AddSecondary(G4Track *aSecondary)
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double GetPDGCharge() const
G4ThreeVector GetDeltaPosition() const
G4double GetScintillationYieldByParticleType(const G4Track &aTrack, const G4Step &aStep)
#define ns
Definition: xmlparse.cc:614
G4MaterialPropertyVector * GetProperty(const char *key)
const G4TouchableHandle & GetTouchableHandle() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scintillation::RemoveSaturation ( )
inline

Definition at line 355 of file G4Scintillation.hh.

356 {
357  fEmSaturation = nullptr;
358 }

Here is the caller graph for this function:

void G4Scintillation::SetFiniteRiseTime ( const G4bool  state)
inline

Definition at line 301 of file G4Scintillation.hh.

302 {
303  fFiniteRiseTime = state;
304 }

Here is the caller graph for this function:

void G4Scintillation::SetScintillationByParticleType ( const G4bool  scintType)

Definition at line 654 of file G4Scintillation.cc.

655 {
656  if (fEmSaturation) {
657  G4Exception("G4Scintillation::SetScintillationByParticleType", "Scint02",
658  JustWarning, "Redefinition: Birks Saturation is replaced by ScintillationByParticleType!");
660  }
661  fScintillationByParticleType = scintType;
662 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scintillation::SetScintillationExcitationRatio ( const G4double  ratio)
inline

Definition at line 325 of file G4Scintillation.hh.

326 {
327  fExcitationRatio = ratio;
328 }

Here is the caller graph for this function:

void G4Scintillation::SetScintillationTrackInfo ( const G4bool  trackType)
inline

Definition at line 373 of file G4Scintillation.hh.

374 {
375  fScintillationTrackInfo = trackType;
376 }

Here is the caller graph for this function:

void G4Scintillation::SetScintillationYieldFactor ( const G4double  yieldfactor)
inline

Definition at line 313 of file G4Scintillation.hh.

314 {
315  fYieldFactor = yieldfactor;
316 }

Here is the caller graph for this function:

void G4Scintillation::SetStackPhotons ( const G4bool  stackingFlag)
inline

Definition at line 385 of file G4Scintillation.hh.

386 {
387  fStackingFlag = stackingFlag;
388 }

Here is the caller graph for this function:

void G4Scintillation::SetTrackSecondariesFirst ( const G4bool  state)
inline

Definition at line 289 of file G4Scintillation.hh.

290 {
291  fTrackSecondariesFirst = state;
292 }

Here is the caller graph for this function:

Member Data Documentation

G4PhysicsTable* G4Scintillation::fFastIntegralTable
protected

Definition at line 241 of file G4Scintillation.hh.

G4PhysicsTable* G4Scintillation::fSlowIntegralTable
protected

Definition at line 242 of file G4Scintillation.hh.


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